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