PandaRoot
PndFTSCATrackParam.h
Go to the documentation of this file.
1 //****************************************************************************
2 //* This file is part of PandaRoot. *
3 //* *
4 //* PandaRoot is distributed under the terms of the *
5 //* GNU General Public License (GPL) version 3, *
6 //* copied verbatim in the file "LICENSE". *
7 //* *
8 //* Copyright (C) 2006 - 2024 FAIR GmbH and copyright holders of PandaRoot *
9 //* The copyright holders are listed in the file "COPYRIGHTHOLDERS". *
10 //* The authors are listed in the file "AUTHORS". *
11 //****************************************************************************
12 
13 
14 //-*- Mode: C++ -*-
15 // $Id: PndFTSCATrackParam.h,v 1.5 2011/05/20 16:28:05 fisyak Exp $
16 // ************************************************************************
17 // This file is property of and copyright by the ALICE HLT Project *
18 // ALICE Experiment at CERN, All rights reserved. *
19 // See cxx source for full Copyright notice *
20 // *
21 //*************************************************************************
22 
23 #ifndef PNDFTSCATRACKPARAM_H
24 #define PNDFTSCATRACKPARAM_H
25 #include <string.h>
26 #include "PndFTSCAMath.h"
28 #include "PndFTSCADef.h"
29 #include "FTSCAHits.h"
30 
31 #ifdef PANDA_FTS
32 
40 class PndFTSCATrackParam {
41  friend std::istream &operator>>(std::istream &, PndFTSCATrackParam &);
42  friend std::ostream &operator<<(std::ostream &, const PndFTSCATrackParam &);
43 
44  public:
45  PndFTSCATrackParam() { Reset(); }
47  : x(v.X()[i]), y(v.Y()[i]), tx(v.Tx()[i]), ty(v.Ty()[i]), qp(v.QP()[i]), z(v.Z()[i]), C00(v.Cov(0)[i]), C10(v.Cov(1)[i]), C11(v.Cov(2)[i]), C20(v.Cov(3)[i]), C21(v.Cov(4)[i]),
48  C22(v.Cov(5)[i]), C30(v.Cov(6)[i]), C31(v.Cov(7)[i]), C32(v.Cov(8)[i]), C33(v.Cov(9)[i]), C40(v.Cov(10)[i]), C41(v.Cov(11)[i]), C42(v.Cov(12)[i]), C43(v.Cov(13)[i]),
49  C44(v.Cov(14)[i]), chi2(v.Chi2()[i]), ndf(v.NDF()[i]), fAlpha(v.Angle()[i])
50  {
51  }
52 
53  float X0() const { return z; }
54  float X1() const { return x; }
55  float X2() const { return y; }
56  float Tx1() const { return tx; }
57  float Tx2() const { return ty; }
58 
59  float X() const { return x; }
60  float Y() const { return y; }
61  float Z() const { return z; }
62  float Tx() const { return tx; }
63  float Ty() const { return ty; }
64  float QP() const { return qp; }
65 
66  float Chi2() const { return chi2; }
67  int NDF() const { return ndf; }
68 
69  float QMomentum() const { return QP(); } // used for triplets comparison
70  float Angle() const { return fAlpha; }
71 
72  float Err2X1() const { return Err2X(); }
73  float Err2X2() const { return Err2Y(); }
74  float Err2Tx1() const { return Err2Tx(); }
75  float Err2Tx2() const { return Err2Ty(); }
76 
77  float Err2X() const { return C00; }
78  float Err2Y() const { return C11; }
79  float Err2Tx() const { return C22; }
80  float Err2Ty() const { return C33; }
81  float Err2QP() const { return C44; }
82  float Err2QMomentum() const { return Err2QP(); }
83 
84  float *Par() { return &x; }
85  float *Cov() { return &C00; }
86 
87  float Par(int i) const
88  {
89  switch (i) {
90  case 0: return x;
91  case 1: return y;
92  case 2: return tx;
93  case 3: return ty;
94  case 4: return qp;
95  }
96  assert(0);
97  return -1;
98  }
99 
100  float Cov(int i) const
101  {
102  switch (i) {
103  case 0: return C00;
104  case 1: return C10;
105  case 2: return C11;
106  case 3: return C20;
107  case 4: return C21;
108  case 5: return C22;
109  case 6: return C30;
110  case 7: return C31;
111  case 8: return C32;
112  case 9: return C33;
113  case 10: return C40;
114  case 11: return C41;
115  case 12: return C42;
116  case 13: return C43;
117  case 14: return C44;
118  }
119  assert(0);
120  return -1;
121  }
122 
123  void SetCov(int i, float v)
124  {
125  switch (i) {
126  case 0: C00 = v; break;
127  case 1: C10 = v; break;
128  case 2: C11 = v; break;
129  case 3: C20 = v; break;
130  case 4: C21 = v; break;
131  case 5: C22 = v; break;
132  case 6: C30 = v; break;
133  case 7: C31 = v; break;
134  case 8: C32 = v; break;
135  case 9: C33 = v; break;
136  case 10: C40 = v; break;
137  case 11: C41 = v; break;
138  case 12: C42 = v; break;
139  case 13: C43 = v; break;
140  case 14: C44 = v; break;
141  }
142  }
143 
144  void SetX(float v) { x = v; }
145  void SetY(float v) { y = v; }
146  void SetZ(float v) { z = v; }
147 
148  void SetTX(float v) { tx = v; }
149  void SetTY(float v) { ty = v; }
150  void SetQP(float v) { qp = v; }
151 
152  void SetChi2(float v) { chi2 = v; }
153  void SetNDF(int v) { ndf = v; }
154 
155  void Print() const;
156 
157  bool TransportToX0(const float &x_, const float &Bz)
158  {
159  UNUSED_PARAM1(Bz);
160  return TransportToX0Line(x_);
161  };
162 
163  bool Transport(const FTSCAHit &hit, const PndFTSCAParam &param)
164  {
166  v.ConvertTrackParamToVector(this, 1);
167  const float_m &mask = v.Transport(hit, param);
168  *this = PndFTSCATrackParam(v, 0);
169  return mask[0];
170  };
171 
172  bool Filter(const FTSCAHit &hit, const PndFTSCAParam &param)
173  {
175  v.ConvertTrackParamToVector(this, 1);
176  const float_m &mask = v.Filter(hit, param);
177  *this = PndFTSCATrackParam(v, 0);
178  return mask[0];
179  };
180 
181  void Reset()
182  {
183  x = y = tx = ty = qp = z = C00 = C10 = C11 = C20 = C21 = C22 = C30 = C31 = C32 = C33 = C40 = C41 = C42 = C43 = C44 = 0.f;
184  ndf = 0;
185  chi2 = -1;
186  }
187 
188  bool TransportToX0Line(float x0);
189 
190  bool IsValid() const { return chi2 != -1; }
191  void SetAsInvalid() { chi2 = -1; }
192 
193  bool Rotate(float) { return 1; }; // don't need rotation in CBM
194 
195  private:
196  float x, y, tx, ty, qp, z, C00, C10, C11, C20, C21, C22, C30, C31, C32, C33, C40, C41, C42, C43, C44;
197 
198  float chi2; // the chi^2 value
199  short ndf; // the Number of Degrees of Freedom
200 
201  float fAlpha; // coor system
202 };
203 
204 inline bool PndFTSCATrackParam::TransportToX0Line(float x0_out)
205 {
206  float dz = (x0_out - z);
207 
208  x += tx * dz;
209  y += ty * dz;
210  z += dz;
211 
212  const float dzC32_in = dz * C32;
213 
214  C21 += dzC32_in;
215  C10 += dz * (C21 + C30);
216 
217  const float C20_in = C20;
218 
219  C20 += dz * C22;
220  C00 += dz * (C20 + C20_in);
221 
222  const float C31_in = C31;
223 
224  C31 += dz * C33;
225  C11 += dz * (C31 + C31_in);
226  C30 += dzC32_in;
227 
228  C40 += dz * C42;
229  C41 += dz * C43;
230 
231  return 1;
232 }
233 
234 #else // PANDA_FTS
235 
237  friend std::istream &operator>>(std::istream &, PndFTSCATrackParam &);
238  friend std::ostream &operator<<(std::ostream &, const PndFTSCATrackParam &);
239 
240  public:
241  PndFTSCATrackParam() { Reset(); }
242  PndFTSCATrackParam(const TrackParamVector &v, int i) : fX(v.X()[i]), fSignCosPhi(v.SignCosPhi()[i]), fChi2(v.Chi2()[i]), fNDF(v.NDF()[i]), fAlpha(v.Angle()[i])
243  {
244  for (int j = 0; j < 5; ++j)
245  fP[j] = v.Par()[j][i];
246  for (int j = 0; j < 15; ++j)
247  fC[j] = v.Cov()[j][i];
248  }
249 
250  float X0() const { return fX; }
251  float X1() const { return Y(); }
252  float X2() const { return Z(); }
253 
254  float X() const { return fX; }
255  float Y() const { return fP[0]; }
256  float Z() const { return fP[1]; }
257  float SinPhi() const { return fP[2]; }
258  float DzDs() const { return fP[3]; }
259  float QPt() const { return fP[4]; }
260  float QP() const { return QPt() / sqrt(DzDs() * DzDs() + 1); }
261  float QMomentum() const { return QPt(); }
262 
263  float SignCosPhi() const { return fSignCosPhi; }
264  float Chi2() const { return fChi2; }
265  int NDF() const { return fNDF; }
266 
267  float Err2X1() const { return Err2Y(); }
268  float Err2X2() const { return Err2Z(); }
269 
270  float Err2Y() const { return fC[0]; }
271  float Err2Z() const { return fC[2]; }
272  float Err2SinPhi() const { return fC[5]; }
273  float Err2DzDs() const { return fC[9]; }
274  float Err2QPt() const { return fC[14]; }
275 
276  float Angle() const { return fAlpha; }
277 
278  float Kappa(float Bz) const { return fP[4] * Bz; }
279  float CosPhi() const { return fSignCosPhi * CAMath::Sqrt(1 - SinPhi() * SinPhi()); }
280 
281  float Err2QMomentum() const { return fC[14]; }
282 
283  const float *Par() const { return fP; }
284  const float *Cov() const { return fC; }
285 
286  bool GetXYZPxPyPzQ(float &x, float &y, float &z, float &px, float &py, float &pz, int &q, float cov[21]) const;
287 
288  void SetSinPhi(float v) { fP[2] = v; }
289 
290  void GetDCAPoint(float x, float y, float z, float &px, float &py, float &pz, float Bz) const;
291 
292  bool TransportToX0(float x, float Bz, float maxSinPhi = .999);
293 
294  bool Rotate(float alpha, float maxSinPhi = .999);
295 
296  bool Transport(const FTSCAHit &hit, const PndFTSCAParam &param);
297 
298  bool IsValid() const { return fChi2 != -1; }
299  void SetAsInvalid() { fChi2 = -1; }
300 
301  private:
302  void Reset()
303  {
304  fX = 0;
305  fSignCosPhi = 0;
306  for (int i = 0; i < 5; i++)
307  fP[i] = 0;
308  for (int i = 0; i < 15; i++)
309  fC[i] = 0;
310  fChi2 = 0;
311  fNDF = 0;
312  }
313  float S(float x, float y, float Bz) const;
314 
315  float fX; // x position
316  float fSignCosPhi; // sign of cosPhi // phi = arctg (Dy/Dx)
317  float fP[5]; // 'active' track parameters: Y, Z, SinPhi = dy/sqrt(dx*dx + dy*dy), DzDs = dz/sqrt(dx*dx + dy*dy), q/Pt
318  float fC[15]; // the covariance matrix for Y,Z,SinPhi,..
319  float fChi2; // the chi^2 value
320  int fNDF; // the Number of Degrees of Freedom
321 
322  float fAlpha; // coor system
323 };
324 
325 inline bool PndFTSCATrackParam::TransportToX0(float x, float Bz, float maxSinPhi)
326 {
327  //* Transport the track parameters to X=x, using linearization at t0, and the field value Bz
328  //* maxSinPhi is the max. allowed value for |t0.SinPhi()|
329  //* linearisation of trajectory t0 is also transported to X=x,
330  //* returns 1 if OK
331  //*
332 
333  const float ex = CosPhi();
334  const float ey = SinPhi();
335  const float k = QPt() * Bz;
336  const float dx = x - X();
337 
338  const float ey1 = k * dx + ey;
339 
340  // check for intersection with X=x
341 
342  if (CAMath::Abs(ey1) > maxSinPhi)
343  return 0;
344 
345  float ex1 = CAMath::Sqrt(1.f - ey1 * ey1);
346  if (ex < 0)
347  ex1 = -ex1;
348 
349  const float dx2 = dx * dx;
350  const float ss = ey + ey1;
351  const float cc = ex + ex1;
352 
353  if ((CAMath::Abs(cc) < 1.e-4 || CAMath::Abs(ex) < 1.e-4 || CAMath::Abs(ex1) < 1.e-4))
354  return 0;
355 
356  const float cci = 1.f / cc;
357  const float exi = 1.f / ex;
358  const float ex1i = 1.f / ex1;
359 
360  const float tg = ss * cci; // tan((phi1+phi)/2)
361 
362  const float dy = dx * tg;
363  float dl = dx * CAMath::Sqrt(1.f + tg * tg);
364 
365  if (cc < 0)
366  dl = -dl;
367  float dSin = dl * k * 0.5;
368  if (dSin > 1.f)
369  dSin = 1.f;
370  if (dSin < -1.f)
371  dSin = -1.f;
372  const float dS = (CAMath::Abs(k) > 1.e-4) ? (2 * CAMath::ASin(dSin) / k) : dl;
373  const float dz = dS * DzDs();
374 
375  // float H0[5] = { 1,0, h2, 0, h4 };
376  // float H1[5] = { 0, 1, 0, dS, 0 };
377  // float H2[5] = { 0, 0, 1, 0, dxBz };
378  // float H3[5] = { 0, 0, 0, 1, 0 };
379  // float H4[5] = { 0, 0, 0, 0, 1 };
380 
381  const float h2 = dx * (1.f + ey * ey1 + ex * ex1) * exi * ex1i * cci;
382  const float h4 = dx2 * (cc + ss * ey1 * ex1i) * cci * cci * Bz;
383  const float dxBz = dx * Bz;
384 
385  fX = X() + dx;
386  fP[0] = Y() + dy;
387  fP[1] = Z() + dz;
388  fP[2] = ey1;
389 
390 #if 1
391  const float c00 = fC[0];
392  const float c10 = fC[1];
393  const float c11 = fC[2];
394  const float c20 = fC[3];
395  const float c21 = fC[4];
396  const float c22 = fC[5];
397  const float c30 = fC[6];
398  const float c31 = fC[7];
399  const float c32 = fC[8];
400  const float c33 = fC[9];
401  const float c40 = fC[10];
402  const float c41 = fC[11];
403  const float c42 = fC[12];
404  const float c43 = fC[13];
405  const float c44 = fC[14];
406 
407  fC[0] = (c00 + h2 * h2 * c22 + h4 * h4 * c44 + 2.f * (h2 * c20 + h4 * c40 + h2 * h4 * c42));
408 
409  fC[1] = c10 + h2 * c21 + h4 * c41 + dS * (c30 + h2 * c32 + h4 * c43);
410  fC[2] = c11 + 2.f * dS * c31 + dS * dS * c33;
411 
412  fC[3] = c20 + h2 * c22 + h4 * c42 + dxBz * (c40 + h2 * c42 + h4 * c44);
413  fC[4] = c21 + dS * c32 + dxBz * (c41 + dS * c43);
414  fC[5] = c22 + 2.f * dxBz * c42 + dxBz * dxBz * c44;
415 
416  fC[6] = c30 + h2 * c32 + h4 * c43;
417  fC[7] = c31 + dS * c33;
418  fC[8] = c32 + dxBz * c43;
419  fC[9] = c33;
420 
421  fC[10] = c40 + h2 * c42 + h4 * c44;
422  fC[11] = c41 + dS * c43;
423  fC[12] = c42 + dxBz * c44;
424  fC[13] = c43;
425  fC[14] = c44;
426 #else
427  fC[0] = (fC[0] + h2 * h2 * fC[5] + h4 * h4 * fC[14] + 2 * (h2 * fC[3] + h4 * fC[10] + h2 * h4 * fC[12]));
428 
429  fC[1] = fC[1] + h2 * fC[4] + h4 * fC[11] + dS * (fC[6] + h2 * fC[8] + h4 * fC[13]);
430  fC[2] = fC[2] + 2 * dS * fC[7] + dS * dS * fC[9];
431 
432  fC[3] = fC[3] + h2 * fC[5] + h4 * fC[12] + dxBz * (fC[10] + h2 * fC[12] + h4 * fC[14]);
433  fC[4] = fC[4] + dS * fC[8] + dxBz * (fC[11] + dS * fC[13]);
434  fC[5] = fC[5] + 2 * dxBz * fC[12] + dxBz * dxBz * fC[14];
435 
436  fC[6] = fC[6] + h2 * fC[8] + h4 * fC[13];
437  fC[7] = fC[7] + dS * fC[9];
438  fC[8] = fC[8] + dxBz * fC[13];
439  fC[9] = fC[9];
440 
441  fC[10] = fC[10] + h2 * fC[12] + h4 * fC[14];
442  fC[11] = fC[11] + dS * fC[13];
443  fC[12] = fC[12] + dxBz * fC[14];
444  fC[13] = fC[13];
445  fC[14] = fC[14];
446 #endif
447 
448  return 1;
449 }
450 
451 #endif // PANDA_FTS
452 
454 
455 std::istream &operator>>(std::istream &in, PndFTSCATrackParam &ot);
456 std::ostream &operator<<(std::ostream &out, const PndFTSCATrackParam &ot);
457 
458 #endif
static T ASin(const T &x)
basic_istream< char, char_traits< char > > istream
bool GetXYZPxPyPzQ(float &x, float &y, float &z, float &px, float &py, float &pz, int &q, float cov[21]) const
float Kappa(float Bz) const
bool Transport(const FTSCAHit &hit, const PndFTSCAParam &param)
bool Rotate(float alpha, float maxSinPhi=.999)
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:40
float Err2QMomentum() const
static T Sqrt(const T &x)
Definition: PndCAMath.h:57
friend std::istream & operator>>(std::istream &, PndFTSCATrackParam &)
friend std::ostream & operator<<(std::ostream &, const PndFTSCATrackParam &)
PndFTSCATrackParam(const TrackParamVector &v, int i)
__m128 v
Definition: P4_F32vec4.h:15
unsigned int i
Definition: P4_F32vec4.h:33
void GetDCAPoint(float x, float y, float z, float &px, float &py, float &pz, float Bz) const
static T Abs(const T &x)
Definition: PndCAMath.h:68
float_m Transport(const int_v &ista, const PndFTSCAParam &param, const float_m &mask=float_m(true))
float Err2SinPhi() const
const float * Par() const
const float * Cov() const
PndFTSCATrackParam TrackParam
bool TransportToX0(float x, float Bz, float maxSinPhi=.999)
float f
Definition: P4_F32vec4.h:32
basic_ostream< char, char_traits< char > > ostream
double alpha
Definition: f_Init.h:31
const float_v & Cov(int i) const
void ConvertTrackParamToVector(PndFTSCATrackParam t0[float_v::Size], int nTracksV)
double pz[39]
Definition: pipisigmas.h:25
float SignCosPhi() const
const float_v & Par(int i) const
float_m Filter(const FTSCAHitV &hit, const PndFTSCAParam &param, const float_m &mask=float_m(true), const float_v &chi2Cut=10e10f)