FairRoot/PandaRoot
KFParticleSIMD.cxx
Go to the documentation of this file.
1 //----------------------------------------------------------------------------
2 // Implementation of the KFParticleSIMD class
3 // .
4 // @author S.Gorbunov, I.Kisel, I.Kulakov, M.Zyzak
5 // @version 1.0
6 // @since 13.05.07
7 //
8 // Class to reconstruct and store the decayed particle parameters.
9 // The method is described in CBM-SOFT note 2007-003,
10 // ``Reconstruction of decayed particles based on the Kalman filter'',
11 // http://www.gsi.de/documents/DOC-2007-May-14-1.pdf
12 //
13 // This class is ALICE interface to general mathematics in KFParticleSIMDCore
14 //
15 // -= Copyright &copy ALICE HLT and CBM L1 Groups =-
16 //____________________________________________________________________________
17 
18 
19 #include "KFParticleSIMD.h"
20 #include "KFParticle.h"
21 #include "KFParticleDatabase.h"
22 
23 #ifdef HomogeneousField
24 #include "KFPTrack.h"
25 #include "KFPVertex.h"
26 #endif
27 
28 #ifdef NonhomogeneousField
29 #include "CbmKFTrack.h"
30 #endif
31 
32 #ifdef HomogeneousField
33 fvec KFParticleSIMD::fgBz = -5.; //* Bz compoment of the magnetic field
34 #endif
35 
36 static const fvec Zero = 0;
37 static const fvec One = 1;
38 
40 {
41  if (!gamma) {
42  KFParticleSIMD mother;
43  mother+= d1;
44  mother+= d2;
45  *this = mother;
46  } else
47  ConstructGamma(d1, d2);
48 }
49 
50 void KFParticleSIMD::Create( const fvec Param[], const fvec Cov[], fvec Charge, fvec mass /*Int_t PID*/ )
51 {
52  // Constructor from "cartesian" track, PID hypothesis should be provided
53  //
54  // Param[6] = { X, Y, Z, Px, Py, Pz } - position and momentum
55  // Cov [21] = lower-triangular part of the covariance matrix:
56  //
57  // ( 0 . . . . . )
58  // ( 1 2 . . . . )
59  // Cov. matrix = ( 3 4 5 . . . ) - numbering of covariance elements in Cov[]
60  // ( 6 7 8 9 . . )
61  // ( 10 11 12 13 14 . )
62  // ( 15 16 17 18 19 20 )
63  fvec C[21];
64  for( int i=0; i<21; i++ ) C[i] = Cov[i];
65 
66 // TParticlePDG* particlePDG = TDatabasePDG::Instance()->GetParticle(PID);
67 // fvec mass = (particlePDG) ? particlePDG->Mass() :0.13957;
68 
69  KFParticleBaseSIMD::Initialize( Param, C, Charge, mass );
70 }
71 
72 #ifdef NonhomogeneousField
73 KFParticleSIMD::KFParticleSIMD( CbmKFTrackInterface* Track[], int NTracks, Int_t *qHypo, const Int_t *pdg)
74 {
75  Create(Track, NTracks, qHypo, pdg);
76 }
77 
78 void KFParticleSIMD::Create(CbmKFTrackInterface* Track[], int NTracks, Int_t *qHypo, const Int_t *pdg)
79 {
80  fvec m[6];
81  fvec V[15];
82  fvec TrMass;
83 
84  fDaughterIds.push_back( fvec(-1) );
85  for (int j=0; j<NTracks; j++)
86  {
87  fDaughterIds.back()[j] = Track[j]->Id();
88 
89 // TrMass[j] = pdg ? TDatabasePDG::Instance()->GetParticle(*pdg)->Mass() : Track[j]->GetMass();
90  TrMass[j] = pdg ? KFParticleDatabase::Instance()->GetMass(*pdg) : Track[j]->GetMass();
91 
92  fNDF[j] = Track[j]->GetRefNDF();
93  fChi2[j] = Track[j]->GetRefChi2();
94 
95  for(int i = 0; i < 6; i++) m[i][j] = Track[j]->GetTrack()[i];
96  for(int i = 0; i < 15; i++) V[i][j] = Track[j]->GetCovMatrix()[i];
97 
98  }
99 
100  fvec a = m[2], b = m[3], qp = m[4];
101 
102  //* convert the track to (x,y,px,py,pz,e,t/p) parameterization
103 
104  fvec c2 = 1./(1. + a*a + b*b);
105  fvec pq = 1./qp;
106  fvec p2 = pq*pq;
107  fvec pz = sqrt(p2*c2);
108  fvec px = a*pz;
109  fvec py = b*pz;
110  fvec eE = sqrt( TrMass*TrMass + p2 );
111 
112  fvec H[3] = { -px*c2, -py*c2, -pz*pq };
113  fvec HE = -pq*p2/eE;
114 
115  fP[0] = m[0];
116  fP[1] = m[1];
117  fP[2] = m[5];
118  fP[3] = px;
119  fP[4] = py;
120  fP[5] = pz;
121  fP[6] = eE;
122  fP[7] = 0;
123 
124  fvec cxpz = H[0]*V[ 3] + H[1]*V[ 6] + H[2]*V[10];
125  fvec cypz = H[0]*V[ 4] + H[1]*V[ 7] + H[2]*V[11];
126  fvec capz = H[0]*V[ 5] + H[1]*V[ 8] + H[2]*V[12];
127  fvec cbpz = H[0]*V[ 8] + H[1]*V[ 9] + H[2]*V[13];
128  fvec cqpz = H[0]*V[12] + H[1]*V[13] + H[2]*V[14];
129  fvec cpzpz = H[0]*H[0]*V[5] +H[1]*H[1]*V[9] + H[2]*H[2]*V[14]
130  + 2*( H[0]*H[1]*V[8] +H[0]*H[2]*V[12] +H[1]*H[2]*V[13]);
131 
132  fC[ 0] = V[0];
133  fC[ 1] = V[1];
134  fC[ 2] = V[2];
135  fC[ 3] = 0.f;
136  fC[ 4] = 0.f;
137  fC[ 5] = 0.f;
138  fC[ 6] = V[3]*pz + a*cxpz;
139  fC[ 7] = V[4]*pz + a*cypz;
140  fC[ 8] = 0.f;
141  fC[ 9] = V[5]*pz*pz + 2*a*pz*capz + a*a*cpzpz;
142  fC[10] = V[6]*pz+b*cxpz;
143  fC[11] = V[7]*pz+b*cypz;
144  fC[12] = 0.f;
145  fC[13] = V[8]*pz*pz + a*pz*cbpz + b*pz*capz + a*b*cpzpz;
146  fC[14] = V[9]*pz*pz + 2*b*pz*cbpz + b*b*cpzpz;
147  fC[15] = cxpz;
148  fC[16] = cypz;
149  fC[17] = 0.f;
150  fC[18] = capz*pz + a*cpzpz;
151  fC[19] = cbpz*pz + b*cpzpz;
152  fC[20] = cpzpz;
153  fC[21] = HE*V[10];
154  fC[22] = HE*V[11];
155  fC[23] = 0.f;
156  fC[24] = HE*( V[12]*pz + a*cqpz );
157  fC[25] = HE*( V[13]*pz + b*cqpz );
158  fC[26] = HE*cqpz;
159  fC[27] = HE*HE*V[14];
160  fC[28] = fC[29] = fC[30] = fC[31] = fC[32] = fC[33] = fC[34] = 0;
161  fC[35] = 1.f;
162 
163  for(int j=0; j<fvecLen; j++)
164  fQ[j] = (qp[j]>0.) ?1 :( (qp[j]<0) ?-1 :0);
166  fIsVtxGuess = false;
167  fIsVtxErrGuess = false;
168  fIsLinearized = false;
169 }
170 
171 KFParticleSIMD::KFParticleSIMD(CbmKFTrackInterface &Track, Int_t *qHypo, const Int_t *pdg)
172 {
173  fDaughterIds.push_back( Track.Id() );
174 
175  fvec m[6];
176  fvec V[15];
177  fvec TrMass;
178 
179  CbmKFTrack Tr(Track);
180  for(unsigned int i=0; i<6; i++)
181  m[i] = Tr.GetTrack()[i];
182  for(unsigned int i=0; i<15; i++)
183  V[i] = Tr.GetCovMatrix()[i];
184 // TrMass = pdg ? TDatabasePDG::Instance()->GetParticle(*pdg)->Mass() :Tr.GetMass();
185  TrMass = pdg ? KFParticleDatabase::Instance()->GetMass(*pdg) :Tr.GetMass();
186 
187  fNDF = Tr.GetRefNDF();
188  fChi2 = Tr.GetRefChi2();
189 
190  fvec a = m[2], b = m[3], qp = m[4];
191 
192  //* convert the track to (x,y,px,py,pz,e,t/p) parameterization
193 
194  fvec c2 = 1./(1. + a*a + b*b);
195  fvec pq = 1./qp;
196  fvec p2 = pq*pq;
197  fvec pz = sqrt(p2*c2);
198  fvec px = a*pz;
199  fvec py = b*pz;
200  fvec eE = sqrt( TrMass*TrMass + p2 );
201 
202  fvec H[3] = { -px*c2, -py*c2, -pz*pq };
203  fvec HE = -pq*p2/eE;
204 
205  fP[0] = m[0];
206  fP[1] = m[1];
207  fP[2] = m[5];
208  fP[3] = px;
209  fP[4] = py;
210  fP[5] = pz;
211  fP[6] = eE;
212  fP[7] = 0;
213 
214  fvec cxpz = H[0]*V[ 3] + H[1]*V[ 6] + H[2]*V[10];
215  fvec cypz = H[0]*V[ 4] + H[1]*V[ 7] + H[2]*V[11];
216  fvec capz = H[0]*V[ 5] + H[1]*V[ 8] + H[2]*V[12];
217  fvec cbpz = H[0]*V[ 8] + H[1]*V[ 9] + H[2]*V[13];
218  fvec cqpz = H[0]*V[12] + H[1]*V[13] + H[2]*V[14];
219  fvec cpzpz = H[0]*H[0]*V[5] +H[1]*H[1]*V[9] + H[2]*H[2]*V[14]
220  + 2*( H[0]*H[1]*V[8] +H[0]*H[2]*V[12] +H[1]*H[2]*V[13]);
221 
222  fC[ 0] = V[0];
223  fC[ 1] = V[1];
224  fC[ 2] = V[2];
225  fC[ 3] = 0.f;
226  fC[ 4] = 0.f;
227  fC[ 5] = 0.f;
228  fC[ 6] = V[3]*pz + a*cxpz;
229  fC[ 7] = V[4]*pz + a*cypz;
230  fC[ 8] = 0.f;
231  fC[ 9] = V[5]*pz*pz + 2*a*pz*capz + a*a*cpzpz;
232  fC[10] = V[6]*pz+b*cxpz;
233  fC[11] = V[7]*pz+b*cypz;
234  fC[12] = 0.f;
235  fC[13] = V[8]*pz*pz + a*pz*cbpz + b*pz*capz + a*b*cpzpz;
236  fC[14] = V[9]*pz*pz + 2*b*pz*cbpz + b*b*cpzpz;
237  fC[15] = cxpz;
238  fC[16] = cypz;
239  fC[17] = 0.f;
240  fC[18] = capz*pz + a*cpzpz;
241  fC[19] = cbpz*pz + b*cpzpz;
242  fC[20] = cpzpz;
243  fC[21] = HE*V[10];
244  fC[22] = HE*V[11];
245  fC[23] = 0.f;
246  fC[24] = HE*( V[12]*pz + a*cqpz );
247  fC[25] = HE*( V[13]*pz + b*cqpz );
248  fC[26] = HE*cqpz;
249  fC[27] = HE*HE*V[14];
250  fC[28] = fC[29] = fC[30] = fC[31] = fC[32] = fC[33] = fC[34] = 0;
251  fC[35] = 1.f;
252 
253  for(int j=0; j<fvecLen; j++)
254  fQ[j] = (qp[j]>0.) ?1 :( (qp[j]<0) ?-1 :0);
256  fIsVtxGuess = false;
257  fIsVtxErrGuess = false;
258  fIsLinearized = false;
259 }
260 
261 // KFParticleSIMD::KFParticleSIMD( CbmKFVertexInterface &vertex )
262 // {
263 // // Constructor from CBM KF vertex
264 // fP[0] = vertex.GetRefX();
265 // fP[1] = vertex.GetRefY();
266 // fP[2] = vertex.GetRefZ();
267 // for( Int_t i=0; i<6; i++)
268 // fC[i] = vertex.GetCovMatrix( )[i];
269 // fChi2 = vertex.GetRefChi2();
270 // fNDF = vertex.GetRefNDF();//2*vertex.GetNContributors() - 3;
271 // fQ = 0;
272 // fAtProductionVertex = 0;
273 // fIsLinearized = 0;
274 // fSFromDecay = 0;
275 // }
276 
277 #endif
278 
279 #ifdef HomogeneousField
281 {
282  // Constructor from ALICE track, PID hypothesis should be provided
283 
284  Double_t r[3];
285  Double_t C[21];
286 
287  for(Int_t iPart = 0; iPart<fvecLen; iPart++)
288  {
289  track[iPart].XvYvZv(r);
290  for(Int_t i=0; i<3; i++)
291  fP[i][iPart] = r[i];
292  track[iPart].PxPyPz(r);
293  for(Int_t i=0; i<3; i++)
294  fP[i+3][iPart] = r[i];
295  fQ[iPart] = track[iPart].Charge();
296  track[iPart].GetCovarianceXYZPxPyPz( C );
297  for(Int_t i=0; i<21; i++)
298  fC[i][iPart] = C[i];
299  }
300 
302  Create(fP,fC,fQ,mass);
303 
304  for(Int_t iPart = 0; iPart<fvecLen; iPart++)
305  {
306  fChi2[iPart] = track[iPart].GetChi2();
307  fNDF[iPart] = track[iPart].GetNDF();
308  }
309 }
310 
311 KFParticleSIMD::KFParticleSIMD(KFPTrack &Track, Int_t *qHypo, const Int_t *pdg)
312 {
313  (void)qHypo; // unused
314  Double_t r[3];
315  Double_t C[21];
316 
317  Track.XvYvZv(r);
318  for(Int_t i=0; i<3; i++)
319  fP[i] = r[i];
320  Track.PxPyPz(r);
321  for(Int_t i=0; i<3; i++)
322  fP[i+3] = r[i];
323  fQ = Track.Charge();
324  Track.GetCovarianceXYZPxPyPz( C );
325  for(Int_t i=0; i<21; i++)
326  fC[i] = C[i];
327 
328  fvec mass = KFParticleDatabase::Instance()->GetMass(*pdg);
329  Create(fP,fC,fQ,mass);
330 
331  fChi2 = Track.GetChi2();
332  fNDF = Track.GetNDF();
333 }
334 
335 KFParticleSIMD::KFParticleSIMD(KFPTrack* Track[], int NTracks, Int_t *qHypo, const Int_t *pdg)
336 {
337  Create(Track, NTracks, qHypo, pdg);
338 }
339 
340 void KFParticleSIMD::Create(KFPTrack* Track[], int NTracks, Int_t *qHypo, const Int_t *pdg)
341 {
342  (void)qHypo; // unused
343  Double_t r[3];
344  Double_t C[21];
345 
346  for(Int_t iPart = 0; iPart<NTracks; iPart++)
347  {
348  Track[iPart]->XvYvZv(r);
349  for(Int_t i=0; i<3; i++)
350  fP[i][iPart] = r[i];
351  Track[iPart]->PxPyPz(r);
352  for(Int_t i=0; i<3; i++)
353  fP[i+3][iPart] = r[i];
354  fQ[iPart] = Track[iPart]->Charge();
355  Track[iPart]->GetCovarianceXYZPxPyPz( C );
356  for(Int_t i=0; i<21; i++)
357  fC[i][iPart] = C[i];
358  }
359 
360  fvec mass = KFParticleDatabase::Instance()->GetMass(*pdg);
361  Create(fP,fC,fQ,mass);
362 
363  for(Int_t iPart = 0; iPart<NTracks; iPart++)
364  {
365  fChi2[iPart] = Track[iPart]->GetChi2();
366  fNDF[iPart] = Track[iPart]->GetNDF();
367  }
368 }
369 
371 {
372  // Constructor from ALICE vertex
373 
374  Double_t r[3];
375  Double_t C[21];
376 
377  vertex.GetXYZ( r );
378  for(Int_t i=0; i<3; i++)
379  fP[i] = r[i];
380  vertex.GetCovarianceMatrix( C );
381  for(Int_t i=0; i<21; i++)
382  fC[i] = C[i];
383  fChi2 = vertex.GetChi2();
384  fNDF = 2*vertex.GetNContributors() - 3;
385  fQ = Zero;
387  fIsLinearized = 0;
388  fSFromDecay = 0;
389 }
390 
391 void KFParticleSIMD::GetExternalTrackParam( const KFParticleBaseSIMD &p, Double_t X[fvecLen], Double_t Alpha[fvecLen], Double_t P[5][fvecLen] )
392 {
393  // Conversion to AliExternalTrackParam parameterization
394 
395  fvec cosA = p.GetPx(), sinA = p.GetPy();
396  fvec pt = sqrt(cosA*cosA + sinA*sinA);
397  fvec pti = Zero;
398  fvec mask = fvec(pt < fvec(1.e-4));
399  pti = (!mask) & (One/pt);
400  cosA = if3(mask, One, cosA*pti);
401  sinA = if3(mask, One, sinA*pti);
402 
403  fvec AlphaSIMD = atan2(sinA,cosA);
404  fvec XSIMD = p.GetX()*cosA + p.GetY()*sinA;
405  fvec PSIMD[5];
406  PSIMD[0]= p.GetY()*cosA - p.GetX()*sinA;
407  PSIMD[1]= p.GetZ();
408  PSIMD[2]= Zero;
409  PSIMD[3]= p.GetPz()*pti;
410  PSIMD[4]= p.GetQ()*pti;
411 
412  for(int iPart=0; iPart<fvecLen; iPart++)
413  {
414  X[iPart] = XSIMD[iPart];
415  Alpha[iPart] = AlphaSIMD[iPart];
416  P[0][iPart] = PSIMD[0][iPart];
417  P[1][iPart] = PSIMD[1][iPart];
418  P[2][iPart] = PSIMD[2][iPart];
419  P[3][iPart] = PSIMD[3][iPart];
420  P[4][iPart] = PSIMD[4][iPart];
421  }
422 }
423 
424 #endif
425 
426 KFParticleSIMD::KFParticleSIMD(KFParticle* parts[], const int nPart)
427 {
428  { // check
429  bool ok = 1;
430  const int nD = (parts[0])->NDaughters();
431  for ( int ie = 1; ie < nPart; ie++ ) {
432  const KFParticle &part = *(parts[ie]);
433  ok &= part.NDaughters() == nD;
434  }
435 // assert(ok);
436  if (!ok) {
437  std::cout << " void CbmKFParticle_simd::Create(CbmKFParticle *parts[], int N) " << std::endl;
438  exit(1);
439  }
440  }
441 
442  fDaughterIds.resize( (parts[0])->NDaughters(), fvec(-1) );
443 
444  for ( int ie = 0; ie < nPart; ie++ ) {
445  KFParticle &part = *(parts[ie]);
446 
447  fId[ie] = part.Id();
448  fDaughterIds.back()[ie] = part.Id();
449 
450  fPDG = part.GetPDG();
451 
452  for( int i = 0; i < 8; ++i )
453  fP[i][ie] = part.Parameters()[i];
454  for( int i = 0; i < 36; ++i )
455  fC[i][ie] = part.CovarianceMatrix()[i];
456 
457  fNDF = part.GetNDF();
458  fChi2[ie] = part.GetChi2();
459  fQ[ie] = part.GetQ();
460  fAtProductionVertex = part.GetAtProductionVertex(); // CHECKME
461 #ifdef NonhomogeneousField
462  L1FieldRegion field(part.GetFieldCoeff());
463  fField.SetOneEntry( ie, field, 0 );
464 #endif
465  }
466 }
467 
469 {
470 
471  fId = part.Id();
472  fNDF = part.GetNDF();
473  fChi2 = part.GetChi2();
474  fQ = part.GetQ();
475  fPDG = part.GetPDG();
477  fIsVtxGuess = 0;
478  fIsVtxErrGuess = 0;
479 
480  SetNDaughters(part.NDaughters());
481  for( int i = 0; i < part.NDaughters(); ++i ) {
482  fDaughterIds.push_back( part.DaughterIds()[i] );
483  }
484 
485  for( int i = 0; i < 8; ++i )
486  fP[i] = part.Parameters()[i];
487  for( int i = 0; i < 36; ++i )
488  fC[i] = part.CovarianceMatrix()[i];
489 
490 #ifdef NonhomogeneousField
491  fField = L1FieldRegion(part.GetFieldCoeff());
492 #endif
493 }
494 
495 fvec KFParticleSIMD::GetDistanceFromVertexXY( const fvec vtx[], const fvec Cv[], fvec &val, fvec &err ) const
496 {
497  //* Calculate DCA distance from vertex (transverse impact parameter) in XY
498  //* v = [xy], Cv=[Cxx,Cxy,Cyy ]-covariance matrix
499 
500  fvec ret = Zero;
501 
502  fvec mP[8];
503  fvec mC[36];
504 
505  Transport( GetDStoPoint(vtx), mP, mC );
506 
507  fvec dx = mP[0] - vtx[0];
508  fvec dy = mP[1] - vtx[1];
509  fvec px = mP[3];
510  fvec py = mP[4];
511  fvec pt = sqrt(px*px + py*py);
512  fvec ex=Zero, ey=Zero;
513  fvec mask = fvec( pt < fvec(1.e-4) );
514  ret = mask;
515  pt = if3(mask, One, pt);
516  ex = (!mask) & (px/pt);
517  ey = (!mask) & (py/pt);
518  val = if3(mask, fvec(1.e4), dy*ex - dx*ey);
519 
520  fvec h0 = -ey;
521  fvec h1 = ex;
522  fvec h3 = (dy*ey + dx*ex)*ey/pt;
523  fvec h4 = -(dy*ey + dx*ex)*ex/pt;
524 
525  err =
526  h0*(h0*GetCovariance(0,0) + h1*GetCovariance(0,1) + h3*GetCovariance(0,3) + h4*GetCovariance(0,4) ) +
527  h1*(h0*GetCovariance(1,0) + h1*GetCovariance(1,1) + h3*GetCovariance(1,3) + h4*GetCovariance(1,4) ) +
528  h3*(h0*GetCovariance(3,0) + h1*GetCovariance(3,1) + h3*GetCovariance(3,3) + h4*GetCovariance(3,4) ) +
529  h4*(h0*GetCovariance(4,0) + h1*GetCovariance(4,1) + h3*GetCovariance(4,3) + h4*GetCovariance(4,4) );
530 
531  if( Cv ){
532  err+= h0*(h0*Cv[0] + h1*Cv[1] ) + h1*(h0*Cv[1] + h1*Cv[2] );
533  }
534 
535  err = sqrt(fabs(err));
536 
537  return ret;
538 }
539 
541 {
542  return GetDistanceFromVertexXY( vtx, 0, val, err );
543 }
544 
545 
547 {
548  //* Calculate distance from vertex [cm] in XY-plane
549 
550  return GetDistanceFromVertexXY( Vtx.fP, Vtx.fC, val, err );
551 }
552 
553 #ifdef HomogeneousField
555 {
556  //* Calculate distance from vertex [cm] in XY-plane
557 
558  return GetDistanceFromVertexXY( KFParticleSIMD(Vtx), val, err );
559 }
560 #endif
561 
563 {
564  //* Calculate distance from vertex [cm] in XY-plane
565  fvec val, err;
566  GetDistanceFromVertexXY( vtx, 0, val, err );
567  return val;
568 }
569 
571 {
572  //* Calculate distance from vertex [cm] in XY-plane
573 
574  return GetDistanceFromVertexXY( Vtx.fP );
575 }
576 
577 #ifdef HomogeneousField
579 {
580  //* Calculate distance from vertex [cm] in XY-plane
581 
583 }
584 #endif
585 
587 {
588  //* Calculate distance to other particle [cm]
589 
590  fvec dS, dS1;
591  GetDStoParticleXY( p, dS, dS1 );
592  fvec mP[8], mC[36], mP1[8], mC1[36];
593  Transport( dS, mP, mC );
594  p.Transport( dS1, mP1, mC1 );
595  fvec dx = mP[0]-mP1[0];
596  fvec dy = mP[1]-mP1[1];
597  return sqrt(dx*dx+dy*dy);
598 }
599 
601 {
602  //* Calculate sqrt(Chi2/ndf) deviation from other particle
603 
604  fvec dS, dS1;
605  GetDStoParticleXY( p, dS, dS1 );
606  fvec mP1[8], mC1[36];
607  p.Transport( dS1, mP1, mC1 );
608 
609  fvec d[2]={ fP[0]-mP1[0], fP[1]-mP1[1] };
610 
611  fvec sigmaS = .1+10.*sqrt( (d[0]*d[0]+d[1]*d[1] )/
612  (mP1[3]*mP1[3]+mP1[4]*mP1[4] ) );
613 
614  fvec h[2] = { mP1[3]*sigmaS, mP1[4]*sigmaS };
615 
616  mC1[0] +=h[0]*h[0];
617  mC1[1] +=h[1]*h[0];
618  mC1[2] +=h[1]*h[1];
619 
620  return GetDeviationFromVertexXY( mP1, mC1 )*sqrt(2./1.);
621 }
622 
623 
624 fvec KFParticleSIMD::GetDeviationFromVertexXY( const fvec vtx[], const fvec Cv[] ) const
625 {
626  //* Calculate sqrt(Chi2/ndf) deviation from vertex
627  //* v = [xyz], Cv=[Cxx,Cxy,Cyy,Cxz,Cyz,Czz]-covariance matrix
628 
629  fvec val, err;
630  fvec problem = GetDistanceFromVertexXY( vtx, Cv, val, err );
631  fvec ret = fvec(1.e4);
632  fvec mask = fvec(problem | fvec(err<fvec(1.e-20)));
633  ret = if3(mask, ret, val/err);
634  return ret;
635 }
636 
637 
639 {
640  //* Calculate sqrt(Chi2/ndf) deviation from vertex
641  //* v = [xyz], Cv=[Cxx,Cxy,Cyy,Cxz,Cyz,Czz]-covariance matrix
642 
643  return GetDeviationFromVertexXY( Vtx.fP, Vtx.fC );
644 }
645 
646 #ifdef HomogeneousField
648 {
649  //* Calculate sqrt(Chi2/ndf) deviation from vertex
650  //* v = [xyz], Cv=[Cxx,Cxy,Cyy,Cxz,Cyz,Czz]-covariance matrix
651 
652  KFParticleSIMD v(Vtx);
653  return GetDeviationFromVertexXY( v.fP, v.fC );
654 }
655 #endif
656 
658 {
659  //* Calculate the opening angle between two particles
660 
661  fvec dS, dS1;
662  GetDStoParticle( p, dS, dS1 );
663  fvec mP[8], mC[36], mP1[8], mC1[36];
664  Transport( dS, mP, mC );
665  p.Transport( dS1, mP1, mC1 );
666  fvec n = sqrt( mP[3]*mP[3] + mP[4]*mP[4] + mP[5]*mP[5] );
667  fvec n1= sqrt( mP1[3]*mP1[3] + mP1[4]*mP1[4] + mP1[5]*mP1[5] );
668  n*=n1;
669  fvec a = Zero;
670  fvec mask = fvec(n>fvec(1.e-8));
671  a = ( mP[3]*mP1[3] + mP[4]*mP1[4] + mP[5]*mP1[5] )/n;
672  a = mask & a;
673  mask = fvec( fabs(a) < One);
674  fvec aPos = fvec(a>=Zero);
675  a = if3(mask, acos(a), if3(aPos, Zero, fvec(3.1415926535)) );
676  return a;
677 }
678 
680 {
681  //* Calculate the opening angle between two particles in XY plane
682 
683  fvec dS, dS1;
684  GetDStoParticleXY( p, dS, dS1 );
685  fvec mP[8], mC[36], mP1[8], mC1[36];
686  Transport( dS, mP, mC );
687  p.Transport( dS1, mP1, mC1 );
688  fvec n = sqrt( mP[3]*mP[3] + mP[4]*mP[4] );
689  fvec n1= sqrt( mP1[3]*mP1[3] + mP1[4]*mP1[4] );
690  n*=n1;
691  fvec a = Zero;
692  fvec mask = fvec(n>fvec(1.e-8));
693  a = ( mP[3]*mP1[3] + mP[4]*mP1[4] )/n;
694  a = mask & a;
695  mask = fvec( fabs(a) < One);
696  fvec aPos = fvec(a>=Zero);
697  a = if3(mask, acos(a), if3(aPos, Zero, fvec(3.1415926535)) );
698  return a;
699 }
700 
702 {
703  //* Calculate the opening angle between two particles in RZ plane
704 
705  fvec dS, dS1;
706  GetDStoParticle( p, dS, dS1 );
707  fvec mP[8], mC[36], mP1[8], mC1[36];
708  Transport( dS, mP, mC );
709  p.Transport( dS1, mP1, mC1 );
710  fvec nr = sqrt( mP[3]*mP[3] + mP[4]*mP[4] );
711  fvec n1r= sqrt( mP1[3]*mP1[3] + mP1[4]*mP1[4] );
712  fvec n = sqrt( nr*nr + mP[5]*mP[5] );
713  fvec n1= sqrt( n1r*n1r + mP1[5]*mP1[5] );
714  n*=n1;
715  fvec a = Zero;
716  fvec mask = fvec(n>fvec(1.e-8));
717  a = ( nr*n1r +mP[5]*mP1[5])/n;
718  a = mask & a;
719  mask = fvec( fabs(a) < One);
720  fvec aPos = fvec(a>=Zero);
721  a = if3(mask, acos(a), if3(aPos, Zero, fvec(3.1415926535)) );
722  return a;
723 }
724 
725 
726 /*
727 
728 #include "AliExternalTrackParam.h"
729 
730 void KFParticleSIMD::GetDStoParticleALICE( const KFParticleSIMDBaseSIMD &p,
731  fvec &DS, fvec &DS1 )
732  const
733 {
734  DS = DS1 = 0;
735  fvec x1, a1, x2, a2;
736  fvec par1[5], par2[5], cov[15];
737  for(int i=0; i<15; i++) cov[i] = 0;
738  cov[0] = cov[2] = cov[5] = cov[9] = cov[14] = .001;
739 
740  GetExternalTrackParam( *this, x1, a1, par1 );
741  GetExternalTrackParam( p, x2, a2, par2 );
742 
743  AliExternalTrackParam t1(x1,a1, par1, cov);
744  AliExternalTrackParam t2(x2,a2, par2, cov);
745 
746  fvec xe1=0, xe2=0;
747  t1.GetDCA( &t2, -GetFieldAlice(), xe1, xe2 );
748  t1.PropagateTo( xe1, -GetFieldAlice() );
749  t2.PropagateTo( xe2, -GetFieldAlice() );
750 
751  fvec xyz1[3], xyz2[3];
752  t1.GetXYZ( xyz1 );
753  t2.GetXYZ( xyz2 );
754 
755  DS = GetDStoPoint( xyz1 );
756  DS1 = p.GetDStoPoint( xyz2 );
757 
758  return;
759 }
760 */
761 
762  // * Pseudo Proper Time of decay = (r*pt) / |pt| * M/|pt|
763 fvec KFParticleSIMD::GetPseudoProperDecayTime( const KFParticleSIMD &pV, const fvec& mass, fvec* timeErr2 ) const
764 { // TODO optimize with respect to time and stability
765  const fvec ipt2 = 1/( Px()*Px() + Py()*Py() );
766  const fvec mipt2 = mass*ipt2;
767  const fvec dx = X() - pV.X();
768  const fvec dy = Y() - pV.Y();
769 
770  if ( timeErr2 ) {
771  // -- calculate error = sigma(f(r)) = f'Cf'
772  // r = {x,y,px,py,x_pV,y_pV}
773  // df/dr = { px*m/pt^2,
774  // py*m/pt^2,
775  // ( x - x_pV )*m*(1/pt^2 - 2(px/pt^2)^2),
776  // ( y - y_pV )*m*(1/pt^2 - 2(py/pt^2)^2),
777  // -px*m/pt^2,
778  // -py*m/pt^2 }
779  const fvec f0 = Px()*mipt2;
780  const fvec f1 = Py()*mipt2;
781  const fvec mipt2derivative = mipt2*(1-2*Px()*Px()*ipt2);
782  const fvec f2 = dx*mipt2derivative;
783  const fvec f3 = -dy*mipt2derivative;
784  const fvec f4 = -f0;
785  const fvec f5 = -f1;
786 
787  const fvec& mC00 = GetCovariance(0,0);
788  const fvec& mC10 = GetCovariance(0,1);
789  const fvec& mC11 = GetCovariance(1,1);
790  const fvec& mC20 = GetCovariance(3,0);
791  const fvec& mC21 = GetCovariance(3,1);
792  const fvec& mC22 = GetCovariance(3,3);
793  const fvec& mC30 = GetCovariance(4,0);
794  const fvec& mC31 = GetCovariance(4,1);
795  const fvec& mC32 = GetCovariance(4,3);
796  const fvec& mC33 = GetCovariance(4,4);
797  const fvec& mC44 = pV.GetCovariance(0,0);
798  const fvec& mC54 = pV.GetCovariance(1,0);
799  const fvec& mC55 = pV.GetCovariance(1,1);
800 
801  *timeErr2 =
802  f5*mC55*f5 +
803  f5*mC54*f4 +
804  f4*mC44*f4 +
805  f3*mC33*f3 +
806  f3*mC32*f2 +
807  f3*mC31*f1 +
808  f3*mC30*f0 +
809  f2*mC22*f2 +
810  f2*mC21*f1 +
811  f2*mC20*f0 +
812  f1*mC11*f1 +
813  f1*mC10*f0 +
814  f0*mC00*f0;
815  }
816  return ( dx*Px() + dy*Py() )*mipt2;
817 }
818 
820 {
821  Part.SetId(static_cast<int>(Id()[iPart]));
822 
823  Part.CleanDaughtersId();
824  Part.SetNDaughters(DaughterIds().size());
825  for( unsigned int i = 0; i < DaughterIds().size(); i++ )
826  Part.AddDaughter(static_cast<int>(DaughterIds()[i][iPart]));
827 
828  Part.SetPDG( static_cast<int>(GetPDG()) );
829 
830  for(int iP=0; iP<8; iP++)
831  Part.Parameters()[iP] = Parameters()[iP][iPart];
832  for(int iC=0; iC<36; iC++)
833  Part.CovarianceMatrix()[iC] = CovarianceMatrix()[iC][iPart];
834 
835  Part.NDF() = static_cast<int>(GetNDF()[iPart]);
836  Part.Chi2() = GetChi2()[iPart];
837  Part.Q() = GetQ()[iPart];
839 #ifdef NonomogeneousField
840  Part.SetFieldCoeff(fField.cx0[iPart], 0);
841  Part.SetFieldCoeff(fField.cx1[iPart], 1);
842  Part.SetFieldCoeff(fField.cx2[iPart], 2);
843  Part.SetFieldCoeff(fField.cy0[iPart], 3);
844  Part.SetFieldCoeff(fField.cy1[iPart], 4);
845  Part.SetFieldCoeff(fField.cy2[iPart], 5);
846  Part.SetFieldCoeff(fField.cz0[iPart], 6);
847  Part.SetFieldCoeff(fField.cz1[iPart], 7);
848  Part.SetFieldCoeff(fField.cz2[iPart], 8);
849  Part.SetFieldCoeff(fField.z0[iPart], 9);
850 #endif
851 }
852 
854 {
855  for(int i=0; i<nPart; i++)
856  GetKFParticle(Part[i],i);
857 }
float GetMass(int pdg)
void GetDStoParticle(const KFParticleSIMD &p, fvec &DS, fvec &DSp) const
friend F32vec4 acos(const F32vec4 &a)
Definition: P4_F32vec4.h:113
PndMultiField * fField
Definition: sim_emc_apd.C:97
std::vector< fvec > fDaughterIds
const fvec & X() const
Double_t p
Definition: anasim.C:58
fvec GetAngleXY(const KFParticleSIMD &p) const
double dy
double r
Definition: RiemannTest.C:14
TObjArray * d
fvec GetChi2() const
Int_t GetNDF() const
Definition: KFParticle.h:497
double mP
Bool_t GetAtProductionVertex() const
Definition: KFParticle.h:121
Int_t i
Definition: run_full.C:25
__m128 m
Definition: P4_F32vec4.h:28
TTree * b
exit(0)
TF1 * f1
Definition: reco_analys2.C:50
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:29
Double_t val[nBoxes][nFEBox]
Definition: createCalib.C:11
const fvec & Y() const
int NDaughters() const
TH1F * h4
fvec GetDeviationFromVertexXY(const fvec v[], const fvec Cv[]=0) const
int GetNContributors() const
Definition: KFPVertex.h:33
Double_t * CovarianceMatrix()
Definition: KFParticle.h:824
int n
c2
Definition: plot_dirc.C:39
TFile * f4
int Charge() const
Definition: KFPTrack.h:62
fvec GetQ() const
Double_t * Parameters()
Definition: KFParticle.h:819
__m128 v
Definition: P4_F32vec4.h:4
fvec GetNDF() const
void GetDStoParticleXY(const KFParticleBaseSIMD &p, fvec &DS, fvec &DSp) const
float GetChi2() const
Definition: KFPVertex.h:31
fvec GetPseudoProperDecayTime(const KFParticleSIMD &primVertex, const fvec &mass, fvec *timeErr2=0) const
fvec GetDStoPoint(const fvec xyz[]) const
TString pt(TString pts, TString exts="px py pz")
Definition: invexp.C:133
int Pic_FED Eff_lEE C()
fvec * CovarianceMatrix()
fvec GetCovariance(int i) const
const Double_t & Chi2() const
Definition: KFParticle.h:138
Int_t a
Definition: anaLmdDigi.C:126
Int_t GetQ() const
Definition: KFParticle.h:487
TFile * f3
FairMCTracks * Track
Definition: drawEveTracks.C:8
void GetXYZ(float *position) const
Definition: KFPVertex.h:17
Double_t
const fvec & Px() const
float_m ok
static KFParticleDatabase * Instance()
PndMCTrack * track
Definition: anaLmdCluster.C:89
fvec GetDeviationFromParticleXY(const KFParticleSIMD &p) const
static const fvec Zero
void SetAtProductionVertex(Bool_t b)
Definition: KFParticle.h:122
F32vec4 fvec
Definition: P4_F32vec4.h:218
void SetNDaughters(int n)
Definition: KFParticle.h:91
std::vector< fvec > & DaughterIds()
friend F32vec4 fabs(const F32vec4 &a)
Definition: P4_F32vec4.h:47
TH1F * h3
const Int_t & Q() const
Definition: KFParticle.h:137
void SetPDG(int pdg)
static const fvec One
void ConstructGamma(const KFParticleSIMD &daughter1, const KFParticleSIMD &daughter2)
TPad * p2
Definition: hist-t7.C:117
double dx
friend F32vec4 atan2(const F32vec4 &y, const F32vec4 &x)
Definition: P4_F32vec4.h:117
double X
Definition: anaLmdDigi.C:68
bool GetCovarianceXYZPxPyPz(float cv[21]) const
Definition: KFPTrack.h:21
int Id() const
const int fvecLen
Definition: P4_F32vec4.h:220
GeV c P
const int & GetPDG() const
const Int_t & NDF() const
Definition: KFParticle.h:139
fvec GetAngle(const KFParticleSIMD &p) const
TFile * f2
Double_t GetChi2() const
Definition: KFParticle.h:492
const std::vector< int > & DaughterIds() const
void PxPyPz(float *position) const
Definition: KFPTrack.h:40
void Transport(fvec dS, fvec P[], fvec C[]) const
static fvec fgBz
int GetPDG() const
void GetKFParticle(KFParticle &Part, int iPart=0)
void GetCovarianceMatrix(float *covmatrix) const
Definition: KFPVertex.h:19
void XvYvZv(float *position) const
Definition: KFPTrack.h:39
void AddDaughter(int id)
Definition: KFParticle.h:92
fvec GetAngleRZ(const KFParticleSIMD &p) const
fvec * Parameters()
int GetNDF() const
Definition: KFPTrack.h:65
fvec GetDistanceFromVertexXY(const fvec vtx[], fvec &val, fvec &err) const
void SetId(int id)
void Create(const fvec Param[], const fvec Cov[], fvec Charge, fvec mass)
#define if3(a, b, c)
Definition: P4_F32vec4.h:91
fvec GetDistanceFromParticleXY(const KFParticleSIMD &p) const
float GetChi2() const
Definition: KFPTrack.h:64
double pz[39]
Definition: pipisigmas.h:14
static void GetExternalTrackParam(const KFParticleBaseSIMD &p, Double_t X[fvecLen], Double_t Alpha[fvecLen], Double_t P[5][fvecLen])
const fvec & Py() const
void CleanDaughtersId()
Definition: KFParticle.h:90