1 //---------------------------------------------------------------------------------
2 // Implementation of the AliKFParticleBase class
4 // @author S.Gorbunov, I.Kisel
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
13 // This class describes general mathematics which is used by AliKFParticle class
15 // -= Copyright © ALICE HLT Group =-
16 //_________________________________________________________________________________
19 #include "AliKFParticleBase.h"
22 ClassImp(AliKFParticleBase)
25 AliKFParticleBase::AliKFParticleBase() :fQ(0), fNDF(-3), fChi2(0), fSFromDecay(0), fAtProductionVertex(0), fIsLinearized(0)
32 void AliKFParticleBase::Initialize( const Double_t Param[], const Double_t Cov[], Int_t Charge, Double_t Mass )
34 // Constructor from "cartesian" track, particle mass hypothesis should be provided
36 // Param[6] = { X, Y, Z, Px, Py, Pz } - position and momentum
37 // Cov [21] = lower-triangular part of the covariance matrix:
41 // Cov. matrix = ( 3 4 5 . . . ) - numbering of covariance elements in Cov[]
43 // ( 10 11 12 13 14 . )
44 // ( 15 16 17 18 19 20 )
47 for( Int_t i=0; i<6 ; i++ ) fP[i] = Param[i];
48 for( Int_t i=0; i<21; i++ ) fC[i] = Cov[i];
50 Double_t energy = TMath::Sqrt( Mass*Mass + fP[3]*fP[3] + fP[4]*fP[4] + fP[5]*fP[5]);
56 fAtProductionVertex = 0;
60 Double_t energyInv = 1./energy;
66 fC[21] = h0*fC[ 6] + h1*fC[10] + h2*fC[15];
67 fC[22] = h0*fC[ 7] + h1*fC[11] + h2*fC[16];
68 fC[23] = h0*fC[ 8] + h1*fC[12] + h2*fC[17];
69 fC[24] = h0*fC[ 9] + h1*fC[13] + h2*fC[18];
70 fC[25] = h0*fC[13] + h1*fC[14] + h2*fC[19];
71 fC[26] = h0*fC[18] + h1*fC[19] + h2*fC[20];
72 fC[27] = ( h0*h0*fC[ 9] + h1*h1*fC[14] + h2*h2*fC[20]
73 + 2*(h0*h1*fC[13] + h0*h2*fC[18] + h1*h2*fC[19] ) );
74 for( Int_t i=28; i<36; i++ ) fC[i] = 0;
78 void AliKFParticleBase::Initialize()
80 //* Initialise covariance matrix and set current parameters to 0.0
82 for( Int_t i=0; i<8; i++) fP[i] = 0;
83 for(Int_t i=0;i<36;++i) fC[i]=0.;
84 fC[0] = fC[2] = fC[5] = 100.;
90 fAtProductionVertex = 0;
91 fVtxGuess[0]=fVtxGuess[1]=fVtxGuess[2]=0.;
95 void AliKFParticleBase::SetVtxGuess( Double_t x, Double_t y, Double_t z )
97 //* Set decay vertex parameters for linearisation
106 Int_t AliKFParticleBase::GetMomentum( Double_t &p, Double_t &error ) const
108 //* Calculate particle momentum
116 Double_t p2 = x2+y2+z2;
118 error = (x2*fC[9]+y2*fC[14]+z2*fC[20] + 2*(x*y*fC[13]+x*z*fC[18]+y*z*fC[19]) );
119 if( error>0 && p>1.e-4 ){
120 error = TMath::Sqrt(error)/p;
126 Int_t AliKFParticleBase::GetPt( Double_t &pt, Double_t &error ) const
128 //* Calculate particle transverse momentum
132 Double_t px2 = px*px;
133 Double_t py2 = py*py;
134 Double_t pt2 = px2+py2;
135 pt = TMath::Sqrt(pt2);
136 error = (px2*fC[9] + py2*fC[14] + 2*px*py*fC[13] );
137 if( error>0 && pt>1.e-4 ){
138 error = TMath::Sqrt(error)/pt;
145 Int_t AliKFParticleBase::GetEta( Double_t &eta, Double_t &error ) const
147 //* Calculate particle pseudorapidity
152 Double_t pt2 = px*px + py*py;
153 Double_t p2 = pt2 + pz*pz;
154 Double_t p = TMath::Sqrt(p2);
160 if( c>1.e-8 ) eta = 0.5*TMath::Log(c);
162 Double_t h3 = -px*pz;
163 Double_t h4 = -py*pz;
164 Double_t pt4 = pt2*pt2;
165 Double_t p2pt4 = p2*pt4;
166 error = (h3*h3*fC[9] + h4*h4*fC[14] + pt4*fC[20] + 2*( h3*(h4*fC[13] + fC[18]*pt2) + pt2*h4*fC[19] ) );
168 if( error>0 && p2pt4>1.e-10 ){
169 error = TMath::Sqrt(error/p2pt4);
177 Int_t AliKFParticleBase::GetPhi( Double_t &phi, Double_t &error ) const
179 //* Calculate particle polar angle
183 Double_t px2 = px*px;
184 Double_t py2 = py*py;
185 Double_t pt2 = px2 + py2;
186 phi = TMath::ATan2(py,px);
187 error = (py2*fC[9] + px2*fC[14] - 2*px*py*fC[13] );
188 if( error>0 && pt2>1.e-4 ){
189 error = TMath::Sqrt(error)/pt2;
196 Int_t AliKFParticleBase::GetR( Double_t &r, Double_t &error ) const
198 //* Calculate distance to the origin
204 r = TMath::Sqrt(x2 + y2);
205 error = (x2*fC[0] + y2*fC[2] - 2*x*y*fC[1] );
206 if( error>0 && r>1.e-4 ){
207 error = TMath::Sqrt(error)/r;
214 Int_t AliKFParticleBase::GetMass( Double_t &m, Double_t &error ) const
216 //* Calculate particle mass
218 // s = sigma^2 of m2/2
220 Double_t s = ( fP[3]*fP[3]*fC[9] + fP[4]*fP[4]*fC[14] + fP[5]*fP[5]*fC[20]
222 +2*( + fP[3]*fP[4]*fC[13] + fP[5]*(fP[3]*fC[18] + fP[4]*fC[19])
223 - fP[6]*( fP[3]*fC[24] + fP[4]*fC[25] + fP[5]*fC[26] ) )
225 Double_t m2 = TMath::Abs(fP[6]*fP[6] - fP[3]*fP[3] - fP[4]*fP[4] - fP[5]*fP[5]);
229 error = TMath::Sqrt(s)/m;
238 Int_t AliKFParticleBase::GetDecayLength( Double_t &l, Double_t &error ) const
240 //* Calculate particle decay length [cm]
249 Double_t p2 = x2+y2+z2;
250 l = t*TMath::Sqrt(p2);
252 error = p2*fC[35] + t*t/p2*(x2*fC[9]+y2*fC[14]+z2*fC[20]
253 + 2*(x*y*fC[13]+x*z*fC[18]+y*z*fC[19]) )
254 + 2*t*(x*fC[31]+y*fC[32]+z*fC[33]);
255 error = TMath::Sqrt(TMath::Abs(error));
262 Int_t AliKFParticleBase::GetLifeTime( Double_t &tauC, Double_t &error ) const
264 //* Calculate particle decay time [s]
268 Double_t cTM = (-fP[3]*fC[31] - fP[4]*fC[32] - fP[5]*fC[33] + fP[6]*fC[34]);
270 error = m*m*fC[35] + 2*fP[7]*cTM + fP[7]*fP[7]*dm*dm;
272 error = TMath::Sqrt( error );
280 void AliKFParticleBase::operator +=( const AliKFParticleBase &Daughter )
282 //* Add daughter via operator+=
284 AddDaughter( Daughter );
287 Double_t AliKFParticleBase::GetSCorrection( const Double_t Part[], const Double_t XYZ[] )
289 //* Get big enough correction for S error to let the particle Part be fitted to XYZ point
291 Double_t d[3] = { XYZ[0]-Part[0], XYZ[1]-Part[1], XYZ[2]-Part[2] };
292 Double_t p2 = Part[3]*Part[3]+Part[4]*Part[4]+Part[5]*Part[5];
293 Double_t sigmaS = (p2>1.e-4) ? ( 10.1+3.*TMath::Sqrt( d[0]*d[0]+d[1]*d[1]+d[2]*d[2]) )/TMath::Sqrt(p2) : 1.;
297 void AliKFParticleBase::GetMeasurement( const Double_t XYZ[], Double_t m[], Double_t V[] ) const
299 //* Get additional covariances V used during measurement
302 GetFieldValue( XYZ, b );
303 const Double_t kCLight = 0.000299792458;
304 b[0]*=kCLight; b[1]*=kCLight; b[2]*=kCLight;
306 Transport( GetDStoPoint(XYZ), m, V );
308 Double_t sigmaS = GetSCorrection( m, XYZ );
315 h[3] = ( h[1]*b[2]-h[2]*b[1] )*GetQ();
316 h[4] = ( h[2]*b[0]-h[0]*b[2] )*GetQ();
317 h[5] = ( h[0]*b[1]-h[1]*b[0] )*GetQ();
346 void AliKFParticleBase::AddDaughter( const AliKFParticleBase &Daughter )
350 if( fNDF<-1 ){ // first daughter -> just copy
352 fQ = Daughter.GetQ();
353 for( Int_t i=0; i<7; i++) fP[i] = Daughter.fP[i];
354 for( Int_t i=0; i<28; i++) fC[i] = Daughter.fC[i];
359 TransportToDecayVertex();
364 if( !fIsLinearized ){
367 GetDStoParticle(Daughter, ds, ds1);
371 Daughter.Transport( ds1, m, mCd );
372 fVtxGuess[0] = .5*( fP[0] + m[0] );
373 fVtxGuess[1] = .5*( fP[1] + m[1] );
374 fVtxGuess[2] = .5*( fP[2] + m[2] );
376 fVtxGuess[0] = fP[0];
377 fVtxGuess[1] = fP[1];
378 fVtxGuess[2] = fP[2];
383 for( Int_t iter=0; iter<maxIter; iter++ ){
386 GetFieldValue( fVtxGuess, b );
387 const Double_t kCLight = 0.000299792458;
388 b[0]*=kCLight; b[1]*=kCLight; b[2]*=kCLight;
391 Double_t *ffP = fP, *ffC = fC, tmpP[8], tmpC[36];
393 GetMeasurement( fVtxGuess, tmpP, tmpC );
398 Double_t m[8], mV[36];
400 if( Daughter.fC[35]>0 ){
401 Daughter.GetMeasurement( fVtxGuess, m, mV );
403 for( Int_t i=0; i<8; i++ ) m[i] = Daughter.fP[i];
404 for( Int_t i=0; i<36; i++ ) mV[i] = Daughter.fC[i];
411 Double_t mSi[6] = { ffC[0]+mV[0],
412 ffC[1]+mV[1], ffC[2]+mV[2],
413 ffC[3]+mV[3], ffC[4]+mV[4], ffC[5]+mV[5] };
415 mS[0] = mSi[2]*mSi[5] - mSi[4]*mSi[4];
416 mS[1] = mSi[3]*mSi[4] - mSi[1]*mSi[5];
417 mS[2] = mSi[0]*mSi[5] - mSi[3]*mSi[3];
418 mS[3] = mSi[1]*mSi[4] - mSi[2]*mSi[3];
419 mS[4] = mSi[1]*mSi[3] - mSi[0]*mSi[4];
420 mS[5] = mSi[0]*mSi[2] - mSi[1]*mSi[1];
422 Double_t s = ( mSi[0]*mS[0] + mSi[1]*mS[1] + mSi[3]*mS[3] );
424 s = ( s > 1.E-20 ) ?1./s :0;
433 //* Residual (measured - estimated)
435 Double_t zeta[3] = { m[0]-ffP[0], m[1]-ffP[1], m[2]-ffP[2] };
440 Double_t mCHt0[7], mCHt1[7], mCHt2[7];
442 mCHt0[0]=ffC[ 0] ; mCHt1[0]=ffC[ 1] ; mCHt2[0]=ffC[ 3] ;
443 mCHt0[1]=ffC[ 1] ; mCHt1[1]=ffC[ 2] ; mCHt2[1]=ffC[ 4] ;
444 mCHt0[2]=ffC[ 3] ; mCHt1[2]=ffC[ 4] ; mCHt2[2]=ffC[ 5] ;
445 mCHt0[3]=ffC[ 6]-mV[ 6]; mCHt1[3]=ffC[ 7]-mV[ 7]; mCHt2[3]=ffC[ 8]-mV[ 8];
446 mCHt0[4]=ffC[10]-mV[10]; mCHt1[4]=ffC[11]-mV[11]; mCHt2[4]=ffC[12]-mV[12];
447 mCHt0[5]=ffC[15]-mV[15]; mCHt1[5]=ffC[16]-mV[16]; mCHt2[5]=ffC[17]-mV[17];
448 mCHt0[6]=ffC[21]-mV[21]; mCHt1[6]=ffC[22]-mV[22]; mCHt2[6]=ffC[23]-mV[23];
450 //* Kalman gain K = mCH'*S
452 Double_t k0[7], k1[7], k2[7];
454 for(Int_t i=0;i<7;++i){
455 k0[i] = mCHt0[i]*mS[0] + mCHt1[i]*mS[1] + mCHt2[i]*mS[3];
456 k1[i] = mCHt0[i]*mS[1] + mCHt1[i]*mS[2] + mCHt2[i]*mS[4];
457 k2[i] = mCHt0[i]*mS[3] + mCHt1[i]*mS[4] + mCHt2[i]*mS[5];
460 //* New estimation of the vertex position
462 if( iter<maxIter-1 ){
463 for(Int_t i=0; i<3; ++i)
464 fVtxGuess[i]= ffP[i] + k0[i]*zeta[0]+k1[i]*zeta[1]+k2[i]*zeta[2];
468 // last itearation -> update the particle
470 //* Add the daughter momentum to the particle momentum
489 //* New estimation of the vertex position r += K*zeta
491 for(Int_t i=0;i<7;++i)
492 fP[i] = ffP[i] + k0[i]*zeta[0] + k1[i]*zeta[1] + k2[i]*zeta[2];
494 //* New covariance matrix C -= K*(mCH')'
496 for(Int_t i=0, k=0;i<7;++i){
497 for(Int_t j=0;j<=i;++j,++k){
498 fC[k] = ffC[k] - (k0[i]*mCHt0[j] + k1[i]*mCHt1[j] + k2[i]*mCHt2[j] );
505 fQ += Daughter.GetQ();
507 fChi2 += (mS[0]*zeta[0] + mS[1]*zeta[1] + mS[3]*zeta[2])*zeta[0]
508 + (mS[1]*zeta[0] + mS[2]*zeta[1] + mS[4]*zeta[2])*zeta[1]
509 + (mS[3]*zeta[0] + mS[4]*zeta[1] + mS[5]*zeta[2])*zeta[2];
515 void AliKFParticleBase::SetProductionVertex( const AliKFParticleBase &Vtx )
517 //* Set production vertex for the particle, when the particle was not used in the vertex fit
519 const Double_t *m = Vtx.fP, *mV = Vtx.fC;
521 Bool_t noS = ( fC[35]<=0 ); // no decay length allowed
524 TransportToDecayVertex();
526 fC[28] = fC[29] = fC[30] = fC[31] = fC[32] = fC[33] = fC[35] = fC[35] = 0;
528 TransportToDS( GetDStoPoint( m ) );
529 fP[7] = -fSFromDecay;
535 InvertSym3( mAi, fC );
539 mB[0][0] = fC[ 6]*mAi[0] + fC[ 7]*mAi[1] + fC[ 8]*mAi[3];
540 mB[0][1] = fC[ 6]*mAi[1] + fC[ 7]*mAi[2] + fC[ 8]*mAi[4];
541 mB[0][2] = fC[ 6]*mAi[3] + fC[ 7]*mAi[4] + fC[ 8]*mAi[5];
543 mB[1][0] = fC[10]*mAi[0] + fC[11]*mAi[1] + fC[12]*mAi[3];
544 mB[1][1] = fC[10]*mAi[1] + fC[11]*mAi[2] + fC[12]*mAi[4];
545 mB[1][2] = fC[10]*mAi[3] + fC[11]*mAi[4] + fC[12]*mAi[5];
547 mB[2][0] = fC[15]*mAi[0] + fC[16]*mAi[1] + fC[17]*mAi[3];
548 mB[2][1] = fC[15]*mAi[1] + fC[16]*mAi[2] + fC[17]*mAi[4];
549 mB[2][2] = fC[15]*mAi[3] + fC[16]*mAi[4] + fC[17]*mAi[5];
551 mB[3][0] = fC[21]*mAi[0] + fC[22]*mAi[1] + fC[23]*mAi[3];
552 mB[3][1] = fC[21]*mAi[1] + fC[22]*mAi[2] + fC[23]*mAi[4];
553 mB[3][2] = fC[21]*mAi[3] + fC[22]*mAi[4] + fC[23]*mAi[5];
555 mB[4][0] = fC[28]*mAi[0] + fC[29]*mAi[1] + fC[30]*mAi[3];
556 mB[4][1] = fC[28]*mAi[1] + fC[29]*mAi[2] + fC[30]*mAi[4];
557 mB[4][2] = fC[28]*mAi[3] + fC[29]*mAi[4] + fC[30]*mAi[5];
559 Double_t z[3] = { m[0]-fP[0], m[1]-fP[1], m[2]-fP[2] };
562 Double_t mAVi[6] = { fC[0]-mV[0], fC[1]-mV[1], fC[2]-mV[2],
563 fC[3]-mV[3], fC[4]-mV[4], fC[5]-mV[5] };
565 if( !InvertSym3( mAVi, mAVi ) ){
567 Double_t dChi2 = ( +(mAVi[0]*z[0] + mAVi[1]*z[1] + mAVi[3]*z[2])*z[0]
568 +(mAVi[1]*z[0] + mAVi[2]*z[1] + mAVi[4]*z[2])*z[1]
569 +(mAVi[3]*z[0] + mAVi[4]*z[1] + mAVi[5]*z[2])*z[2] );
571 // Take Abs(dChi2) here. Negative value of 'det' or 'dChi2' shows that the particle
572 // was not used in the production vertex fit
574 fChi2+= TMath::Abs( dChi2 );
582 fP[3]+= mB[0][0]*z[0] + mB[0][1]*z[1] + mB[0][2]*z[2];
583 fP[4]+= mB[1][0]*z[0] + mB[1][1]*z[1] + mB[1][2]*z[2];
584 fP[5]+= mB[2][0]*z[0] + mB[2][1]*z[1] + mB[2][2]*z[2];
585 fP[6]+= mB[3][0]*z[0] + mB[3][1]*z[1] + mB[3][2]*z[2];
586 fP[7]+= mB[4][0]*z[0] + mB[4][1]*z[1] + mB[4][2]*z[2];
597 d0= mB[0][0]*mV[0] + mB[0][1]*mV[1] + mB[0][2]*mV[3] - fC[ 6];
598 d1= mB[0][0]*mV[1] + mB[0][1]*mV[2] + mB[0][2]*mV[4] - fC[ 7];
599 d2= mB[0][0]*mV[3] + mB[0][1]*mV[4] + mB[0][2]*mV[5] - fC[ 8];
604 fC[ 9]+= d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
606 d0= mB[1][0]*mV[0] + mB[1][1]*mV[1] + mB[1][2]*mV[3] - fC[10];
607 d1= mB[1][0]*mV[1] + mB[1][1]*mV[2] + mB[1][2]*mV[4] - fC[11];
608 d2= mB[1][0]*mV[3] + mB[1][1]*mV[4] + mB[1][2]*mV[5] - fC[12];
613 fC[13]+= d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
614 fC[14]+= d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
616 d0= mB[2][0]*mV[0] + mB[2][1]*mV[1] + mB[2][2]*mV[3] - fC[15];
617 d1= mB[2][0]*mV[1] + mB[2][1]*mV[2] + mB[2][2]*mV[4] - fC[16];
618 d2= mB[2][0]*mV[3] + mB[2][1]*mV[4] + mB[2][2]*mV[5] - fC[17];
623 fC[18]+= d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
624 fC[19]+= d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
625 fC[20]+= d0*mB[2][0] + d1*mB[2][1] + d2*mB[2][2];
627 d0= mB[3][0]*mV[0] + mB[3][1]*mV[1] + mB[3][2]*mV[3] - fC[21];
628 d1= mB[3][0]*mV[1] + mB[3][1]*mV[2] + mB[3][2]*mV[4] - fC[22];
629 d2= mB[3][0]*mV[3] + mB[3][1]*mV[4] + mB[3][2]*mV[5] - fC[23];
634 fC[24]+= d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
635 fC[25]+= d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
636 fC[26]+= d0*mB[2][0] + d1*mB[2][1] + d2*mB[2][2];
637 fC[27]+= d0*mB[3][0] + d1*mB[3][1] + d2*mB[3][2];
639 d0= mB[4][0]*mV[0] + mB[4][1]*mV[1] + mB[4][2]*mV[3] - fC[28];
640 d1= mB[4][0]*mV[1] + mB[4][1]*mV[2] + mB[4][2]*mV[4] - fC[29];
641 d2= mB[4][0]*mV[3] + mB[4][1]*mV[4] + mB[4][2]*mV[5] - fC[30];
646 fC[31]+= d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
647 fC[32]+= d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
648 fC[33]+= d0*mB[2][0] + d1*mB[2][1] + d2*mB[2][2];
649 fC[34]+= d0*mB[3][0] + d1*mB[3][1] + d2*mB[3][2];
650 fC[35]+= d0*mB[4][0] + d1*mB[4][1] + d2*mB[4][2];
654 fC[28] = fC[29] = fC[30] = fC[31] = fC[32] = fC[33] = fC[35] = fC[35] = 0;
656 TransportToDS( fP[7] );
665 void AliKFParticleBase::SetMassConstraint( Double_t Mass, Double_t SigmaMass )
667 //* Set hard( SigmaMass=0 ) or soft (SigmaMass>0) mass constraint
669 Double_t m2 = Mass*Mass; // measurement, weighted by Mass
670 Double_t s2 = m2*SigmaMass*SigmaMass; // sigma^2
672 Double_t p2 = fP[3]*fP[3] + fP[4]*fP[4] + fP[5]*fP[5];
673 Double_t e0 = TMath::Sqrt(m2+p2);
676 mH[0] = mH[1] = mH[2] = 0.;
680 mH[6] = 2*fP[6];//e0;
683 Double_t zeta = e0*e0 - e0*fP[6];
684 zeta = m2 - (fP[6]*fP[6]-p2);
686 Double_t mCHt[8], s2_est=0;
687 for( Int_t i=0; i<8; ++i ){
689 for (Int_t j=0;j<8;++j) mCHt[i] += Cij(i,j)*mH[j];
690 s2_est += mH[i]*mCHt[i];
693 if( s2_est<1.e-20 ) return; // calculated mass error is already 0,
694 // the particle can not be constrained on mass
696 Double_t w2 = 1./( s2 + s2_est );
697 fChi2 += zeta*zeta*w2;
699 for( Int_t i=0, ii=0; i<8; ++i ){
700 Double_t ki = mCHt[i]*w2;
702 for(Int_t j=0;j<=i;++j) fC[ii++] -= ki*mCHt[j];
707 void AliKFParticleBase::SetNoDecayLength()
709 //* Set no decay length for resonances
711 TransportToDecayVertex();
714 h[0] = h[1] = h[2] = h[3] = h[4] = h[5] = h[6] = 0;
717 Double_t zeta = 0 - fP[7];
718 for(Int_t i=0;i<8;++i) zeta -= h[i]*(fP[i]-fP[i]);
723 fChi2 += zeta*zeta*s;
725 for( Int_t i=0, ii=0; i<7; ++i ){
726 Double_t ki = fC[28+i]*s;
728 for(Int_t j=0;j<=i;++j) fC[ii++] -= ki*fC[28+j];
732 fC[28] = fC[29] = fC[30] = fC[31] = fC[32] = fC[33] = fC[35] = fC[35] = 0;
736 void AliKFParticleBase::Construct( const AliKFParticleBase* vDaughters[], Int_t NDaughters,
737 const AliKFParticleBase *Parent, Double_t Mass, Bool_t IsConstrained )
739 //* Full reconstruction in one go
742 bool wasLinearized = fIsLinearized;
743 if( !fIsLinearized || IsConstrained ){
744 //fVtxGuess[0] = fVtxGuess[1] = fVtxGuess[2] = 0; //!!!!
745 fVtxGuess[0] = GetX();
746 fVtxGuess[1] = GetY();
747 fVtxGuess[2] = GetZ();
752 Double_t constraintC[6];
755 for(Int_t i=0;i<6;++i) constraintC[i]=fC[i];
757 for(Int_t i=0;i<6;++i) constraintC[i]=0.;
758 constraintC[0] = constraintC[2] = constraintC[5] = 100.;
762 for( Int_t iter=0; iter<maxIter; iter++ ){
763 fAtProductionVertex = 0;
765 fP[0] = fVtxGuess[0];
766 fP[1] = fVtxGuess[1];
767 fP[2] = fVtxGuess[2];
774 for(Int_t i=0;i<6; ++i) fC[i]=constraintC[i];
775 for(Int_t i=6;i<36;++i) fC[i]=0.;
778 fNDF = IsConstrained ?0 :-3;
782 for( Int_t itr =0; itr<NDaughters; itr++ ){
783 AddDaughter( *vDaughters[itr] );
786 for( Int_t i=0; i<3; i++ ) fVtxGuess[i] = fP[i];
789 fIsLinearized = wasLinearized;
791 if( Mass>=0 ) SetMassConstraint( Mass );
792 if( Parent ) SetProductionVertex( *Parent );
796 void AliKFParticleBase::Convert( bool ToProduction )
798 //* Tricky function - convert the particle error along its trajectory to
799 //* the value which corresponds to its production/decay vertex
800 //* It is done by combination of the error of decay length with the position errors
804 GetFieldValue( fP, fld );
805 const Double_t kCLight = fQ*0.000299792458;
806 fld[0]*=kCLight; fld[1]*=kCLight; fld[2]*=kCLight;
814 if( ToProduction ){ h[0]=-h[0]; h[1]=-h[1]; h[2]=-h[2]; }
815 h[3] = h[1]*fld[2]-h[2]*fld[1];
816 h[4] = h[2]*fld[0]-h[0]*fld[2];
817 h[5] = h[0]*fld[1]-h[1]*fld[0];
821 c = fC[28]+h[0]*fC[35];
822 fC[ 0]+= h[0]*(c+fC[28]);
825 fC[ 1]+= h[1]*fC[28] + h[0]*fC[29];
826 c = fC[29]+h[1]*fC[35];
827 fC[ 2]+= h[1]*(c+fC[29]);
830 fC[ 3]+= h[2]*fC[28] + h[0]*fC[30];
831 fC[ 4]+= h[2]*fC[29] + h[1]*fC[30];
832 c = fC[30]+h[2]*fC[35];
833 fC[ 5]+= h[2]*(c+fC[30]);
836 fC[ 6]+= h[3]*fC[28] + h[0]*fC[31];
837 fC[ 7]+= h[3]*fC[29] + h[1]*fC[31];
838 fC[ 8]+= h[3]*fC[30] + h[2]*fC[31];
839 c = fC[31]+h[3]*fC[35];
840 fC[ 9]+= h[3]*(c+fC[31]);
843 fC[10]+= h[4]*fC[28] + h[0]*fC[32];
844 fC[11]+= h[4]*fC[29] + h[1]*fC[32];
845 fC[12]+= h[4]*fC[30] + h[2]*fC[32];
846 fC[13]+= h[4]*fC[31] + h[3]*fC[32];
847 c = fC[32]+h[4]*fC[35];
848 fC[14]+= h[4]*(c+fC[32]);
851 fC[15]+= h[5]*fC[28] + h[0]*fC[33];
852 fC[16]+= h[5]*fC[29] + h[1]*fC[33];
853 fC[17]+= h[5]*fC[30] + h[2]*fC[33];
854 fC[18]+= h[5]*fC[31] + h[3]*fC[33];
855 fC[19]+= h[5]*fC[32] + h[4]*fC[33];
856 c = fC[33]+h[5]*fC[35];
857 fC[20]+= h[5]*(c+fC[33]);
860 fC[21]+= h[0]*fC[34];
861 fC[22]+= h[1]*fC[34];
862 fC[23]+= h[2]*fC[34];
863 fC[24]+= h[3]*fC[34];
864 fC[25]+= h[4]*fC[34];
865 fC[26]+= h[5]*fC[34];
869 void AliKFParticleBase::TransportToDecayVertex()
871 //* Transport the particle to its decay vertex
873 if( fSFromDecay != 0 ) TransportToDS( -fSFromDecay );
874 if( fAtProductionVertex ) Convert(0);
875 fAtProductionVertex = 0;
878 void AliKFParticleBase::TransportToProductionVertex()
880 //* Transport the particle to its production vertex
882 if( fSFromDecay != -fP[7] ) TransportToDS( -fSFromDecay-fP[7] );
883 if( !fAtProductionVertex ) Convert( 1 );
884 fAtProductionVertex = 1;
888 void AliKFParticleBase::TransportToDS( Double_t dS )
890 //* Transport the particle on dS parameter (SignedPath/Momentum)
892 Transport( dS, fP, fC );
897 Double_t AliKFParticleBase::GetDStoPointLine( const Double_t xyz[] ) const
899 //* Get dS to a certain space point without field
901 Double_t p2 = fP[3]*fP[3] + fP[4]*fP[4] + fP[5]*fP[5];
902 if( p2<1.e-4 ) p2 = 1;
903 return ( fP[3]*(xyz[0]-fP[0]) + fP[4]*(xyz[1]-fP[1]) + fP[5]*(xyz[2]-fP[2]) )/p2;
907 Double_t AliKFParticleBase::GetDStoPointBz( Double_t B, const Double_t xyz[] )
911 //* Get dS to a certain space point for Bz field
912 const Double_t kCLight = 0.000299792458;
913 Double_t bq = B*fQ*kCLight;
914 Double_t pt2 = fP[3]*fP[3] + fP[4]*fP[4];
915 if( pt2<1.e-4 ) return 0;
916 Double_t dx = xyz[0] - fP[0];
917 Double_t dy = xyz[1] - fP[1];
918 Double_t a = dx*fP[3]+dy*fP[4];
921 if( TMath::Abs(bq)<1.e-8 ) dS = a/pt2;
922 else dS = TMath::ATan2( bq*a, pt2 + bq*(dy*fP[3] -dx*fP[4]) )/bq;
929 Double_t ss[2], g[2][5];
933 for( Int_t i=0; i<2; i++){
934 Double_t bs = bq*ss[i];
935 Double_t c = TMath::Cos(bs), s = TMath::Sin(bs);
937 if( TMath::Abs(bq)>1.e-8){
941 const Double_t kOvSqr6 = 1./TMath::Sqrt(6.);
942 sB = (1.-bs*kOvSqr6)*(1.+bs*kOvSqr6)*ss[i];
945 g[i][0] = fP[0] + sB*px + cB*py;
946 g[i][1] = fP[1] - cB*px + sB*py;
947 g[i][2] = fP[2] + ss[i]*pz;
948 g[i][3] = + c*px + s*py;
949 g[i][4] = - s*px + c*py;
954 Double_t dMin = 1.e10;
955 for( Int_t j=0; j<2; j++){
956 Double_t xx = g[j][0]-xyz[0];
957 Double_t yy = g[j][1]-xyz[1];
958 Double_t zz = g[j][2]-xyz[2];
959 Double_t d = xx*xx + yy*yy + zz*zz;
968 Double_t x= g[i][0], y= g[i][1], z= g[i][2], ppx= g[i][3], ppy= g[i][4];
969 Double_t ddx = x-xyz[0];
970 Double_t ddy = y-xyz[1];
971 Double_t ddz = z-xyz[2];
972 Double_t c = ddx*ppx + ddy*ppy + ddz*pz ;
973 Double_t pp2 = ppx*ppx + ppy*ppy + pz*pz;
974 if( TMath::Abs(pp2)>1.e-8 ){
982 void AliKFParticleBase::GetDStoParticleBz( Double_t B, const AliKFParticleBase &p,
983 Double_t &DS, Double_t &DS1 )
986 //* Get dS to another particle for Bz field
991 Double_t px1 = p.fP[3];
992 Double_t py1 = p.fP[4];
993 Double_t pz1 = p.fP[5];
995 const Double_t kCLight = 0.000299792458;
997 Double_t bq = B*fQ*kCLight;
998 Double_t bq1 = B*p.fQ*kCLight;
999 Double_t s=0, ds=0, s1=0, ds1=0;
1001 if( TMath::Abs(bq)>1.e-8 || TMath::Abs(bq1)>1.e-8 ){
1003 Double_t dx = (p.fP[0] - fP[0]);
1004 Double_t dy = (p.fP[1] - fP[1]);
1005 Double_t d2 = (dx*dx+dy*dy);
1007 Double_t p2 = (px *px + py *py);
1008 Double_t p21 = (px1*px1 + py1*py1);
1010 Double_t a = (px*py1 - py*px1);
1011 Double_t b = (px*px1 + py*py1);
1013 Double_t ldx = bq*bq1*dx - bq1*py + bq*py1 ;
1014 Double_t ldy = bq*bq1*dy + bq1*px - bq*px1 ;
1015 Double_t l2 = ldx*ldx + ldy*ldy;
1017 Double_t cS = bq1*p2 + bq*bq1*(dy* px - dx* py) - bq*b;
1018 Double_t cS1= bq*p21 - bq*bq1*(dy*px1 - dx*py1) - bq1*b;
1020 Double_t ca = bq*bq*bq1*d2 +2*( cS + bq*bq*(py1*dx-px1*dy)) ;
1021 Double_t ca1 = bq*bq1*bq1*d2 +2*( cS1 - bq1*bq1*(py*dx-px*dy)) ;
1023 Double_t sa = 4*l2*p2 - ca*ca;
1024 Double_t sa1 = 4*l2*p21 - ca1*ca1;
1029 if( TMath::Abs(bq)>1.e-8){
1030 s = TMath::ATan2( bq*( bq1*(dx*px +dy*py) + a ) , cS )/bq;
1031 ds = TMath::ATan2(TMath::Sqrt(sa),ca)/bq;
1033 s = ( (dx*px + dy*py) + (py*px1-px*py1)/bq1)/p2;
1034 ds = s*s - (d2-2*(px1*dy-py1*dx)/bq1)/p2;
1036 ds = TMath::Sqrt(ds);
1039 if( TMath::Abs(bq1)>1.e-8){
1040 s1 = TMath::ATan2( -bq1*( bq*(dx*px1+dy*py1) + a), cS1 )/bq1;
1041 ds1 = TMath::ATan2(TMath::Sqrt(sa1),ca1)/bq1;
1043 s1 = (-(dx*px1 + dy*py1) + (py*px1-px*py1)/bq)/p21;
1044 ds1 = s1*s1 - (d2+2*(px*dy-py*dx)/bq)/p21;
1045 if( ds1<0 ) ds1 = 0;
1046 ds1 = TMath::Sqrt(ds1);
1050 Double_t ss[2], ss1[2], g[2][5],g1[2][5];
1056 for( Int_t i=0; i<2; i++){
1057 Double_t bs = bq*ss[i];
1058 Double_t c = TMath::Cos(bs), sss = TMath::Sin(bs);
1060 if( TMath::Abs(bq)>1.e-8){
1064 const Double_t kOvSqr6 = 1./TMath::Sqrt(6.);
1065 sB = (1.-bs*kOvSqr6)*(1.+bs*kOvSqr6)*ss[i];
1068 g[i][0] = fP[0] + sB*px + cB*py;
1069 g[i][1] = fP[1] - cB*px + sB*py;
1070 g[i][2] = fP[2] + ss[i]*pz;
1071 g[i][3] = + c*px + sss*py;
1072 g[i][4] = - sss*px + c*py;
1075 c = TMath::Cos(bs); sss = TMath::Sin(bs);
1076 if( TMath::Abs(bq1)>1.e-8){
1080 const Double_t kOvSqr6 = 1./TMath::Sqrt(6.);
1081 sB = (1.-bs*kOvSqr6)*(1.+bs*kOvSqr6)*ss1[i];
1085 g1[i][0] = p.fP[0] + sB*px1 + cB*py1;
1086 g1[i][1] = p.fP[1] - cB*px1 + sB*py1;
1087 g1[i][2] = p.fP[2] + ss[i]*pz1;
1088 g1[i][3] = + c*px1 + sss*py1;
1089 g1[i][4] = - sss*px1 + c*py1;
1094 Double_t dMin = 1.e10;
1095 for( Int_t j=0; j<2; j++){
1096 for( Int_t j1=0; j1<2; j1++){
1097 Double_t xx = g[j][0]-g1[j1][0];
1098 Double_t yy = g[j][1]-g1[j1][1];
1099 Double_t zz = g[j][2]-g1[j1][2];
1100 Double_t d = xx*xx + yy*yy + zz*zz;
1112 Double_t x= g[i][0], y= g[i][1], z= g[i][2], ppx= g[i][3], ppy= g[i][4];
1113 Double_t x1=g1[i1][0], y1= g1[i1][1], z1= g1[i1][2], ppx1= g1[i1][3], ppy1= g1[i1][4];
1117 Double_t a = ppx*ppx1 + ppy*ppy1 + pz*pz1;
1118 Double_t b = dx*ppx1 + dy*ppy1 + dz*pz1;
1119 Double_t c = dx*ppx + dy*ppy + dz*pz ;
1120 Double_t pp2 = ppx*ppx + ppy*ppy + pz*pz;
1121 Double_t pp21= ppx1*ppx1 + ppy1*ppy1 + pz1*pz1;
1122 Double_t det = pp2*pp21 - a*a;
1123 if( TMath::Abs(det)>1.e-8 ){
1124 DS+=(a*b-pp21*c)/det;
1125 DS1+=(a*c-pp2*b)/det;
1132 void AliKFParticleBase::TransportCBM( Double_t dS,
1133 Double_t P[], Double_t C[] ) const
1135 //* Transport the particle on dS, output to P[],C[], for CBM field
1138 TransportLine( dS, P, C );
1142 const Double_t kCLight = 0.000299792458;
1144 Double_t c = fQ*kCLight;
1146 // construct coefficients
1153 Double_t sx=0, sy=0, sz=0, syy=0, syz=0, syyy=0, ssx=0, ssy=0, ssz=0, ssyy=0, ssyz=0, ssyyy=0;
1155 { // get field integrals
1158 Double_t p0[3], p1[3], p2[3];
1160 // line track approximation
1166 p2[0] = fP[0] + px*dS;
1167 p2[1] = fP[1] + py*dS;
1168 p2[2] = fP[2] + pz*dS;
1170 p1[0] = 0.5*(p0[0]+p2[0]);
1171 p1[1] = 0.5*(p0[1]+p2[1]);
1172 p1[2] = 0.5*(p0[2]+p2[2]);
1174 // first order track approximation
1176 GetFieldValue( p0, fld[0] );
1177 GetFieldValue( p1, fld[1] );
1178 GetFieldValue( p2, fld[2] );
1180 Double_t ssy1 = ( 7*fld[0][1] + 6*fld[1][1]-fld[2][1] )*c*dS*dS/96.;
1181 Double_t ssy2 = ( fld[0][1] + 2*fld[1][1] )*c*dS*dS/6.;
1189 GetFieldValue( p0, fld[0] );
1190 GetFieldValue( p1, fld[1] );
1191 GetFieldValue( p2, fld[2] );
1193 sx = c*( fld[0][0] + 4*fld[1][0] + fld[2][0] )*dS/6.;
1194 sy = c*( fld[0][1] + 4*fld[1][1] + fld[2][1] )*dS/6.;
1195 sz = c*( fld[0][2] + 4*fld[1][2] + fld[2][2] )*dS/6.;
1197 ssx = c*( fld[0][0] + 2*fld[1][0])*dS*dS/6.;
1198 ssy = c*( fld[0][1] + 2*fld[1][1])*dS*dS/6.;
1199 ssz = c*( fld[0][2] + 2*fld[1][2])*dS*dS/6.;
1201 Double_t c2[3][3] = { { 5, -4, -1},{ 44, 80, -4},{ 11, 44, 5} }; // /=360.
1202 Double_t cc2[3][3] = { { 38, 8, -4},{ 148, 208, -20},{ 3, 36, 3} }; // /=2520.
1203 for(Int_t n=0; n<3; n++)
1204 for(Int_t m=0; m<3; m++)
1206 syz += c2[n][m]*fld[n][1]*fld[m][2];
1207 ssyz += cc2[n][m]*fld[n][1]*fld[m][2];
1210 syz *= c*c*dS*dS/360.;
1211 ssyz *= c*c*dS*dS*dS/2520.;
1213 syy = c*( fld[0][1] + 4*fld[1][1] + fld[2][1] )*dS;
1214 syyy = syy*syy*syy / 1296;
1217 ssyy = ( fld[0][1]*( 38*fld[0][1] + 156*fld[1][1] - fld[2][1] )+
1218 fld[1][1]*( 208*fld[1][1] +16*fld[2][1] )+
1219 fld[2][1]*( 3*fld[2][1] )
1220 )*dS*dS*dS*c*c/2520.;
1223 fld[0][1]*( fld[0][1]*( 85*fld[0][1] + 526*fld[1][1] - 7*fld[2][1] )+
1224 fld[1][1]*( 1376*fld[1][1] +84*fld[2][1] )+
1225 fld[2][1]*( 19*fld[2][1] ) )+
1226 fld[1][1]*( fld[1][1]*( 1376*fld[1][1] +256*fld[2][1] )+
1227 fld[2][1]*( 62*fld[2][1] ) )+
1228 fld[2][1]*fld[2][1] *( 3*fld[2][1] )
1229 )*dS*dS*dS*dS*c*c*c/90720.;
1234 for( Int_t i=0; i<8; i++ ) for( Int_t j=0; j<8; j++) mJ[i][j]=0;
1236 mJ[0][0]=1; mJ[0][1]=0; mJ[0][2]=0; mJ[0][3]=dS-ssyy; mJ[0][4]=ssx; mJ[0][5]=ssyyy-ssy;
1237 mJ[1][0]=0; mJ[1][1]=1; mJ[1][2]=0; mJ[1][3]=-ssz; mJ[1][4]=dS; mJ[1][5]=ssx+ssyz;
1238 mJ[2][0]=0; mJ[2][1]=0; mJ[2][2]=1; mJ[2][3]=ssy-ssyyy; mJ[2][4]=-ssx; mJ[2][5]=dS-ssyy;
1240 mJ[3][0]=0; mJ[3][1]=0; mJ[3][2]=0; mJ[3][3]=1-syy; mJ[3][4]=sx; mJ[3][5]=syyy-sy;
1241 mJ[4][0]=0; mJ[4][1]=0; mJ[4][2]=0; mJ[4][3]=-sz; mJ[4][4]=1; mJ[4][5]=sx+syz;
1242 mJ[5][0]=0; mJ[5][1]=0; mJ[5][2]=0; mJ[5][3]=sy-syyy; mJ[5][4]=-sx; mJ[5][5]=1-syy;
1243 mJ[6][6] = mJ[7][7] = 1;
1245 P[0] = fP[0] + mJ[0][3]*px + mJ[0][4]*py + mJ[0][5]*pz;
1246 P[1] = fP[1] + mJ[1][3]*px + mJ[1][4]*py + mJ[1][5]*pz;
1247 P[2] = fP[2] + mJ[2][3]*px + mJ[2][4]*py + mJ[2][5]*pz;
1248 P[3] = mJ[3][3]*px + mJ[3][4]*py + mJ[3][5]*pz;
1249 P[4] = mJ[4][3]*px + mJ[4][4]*py + mJ[4][5]*pz;
1250 P[5] = mJ[5][3]*px + mJ[5][4]*py + mJ[5][5]*pz;
1254 MultQSQt( mJ[0], fC, C);
1259 void AliKFParticleBase::TransportBz( Double_t b, Double_t t,
1260 Double_t p[], Double_t e[] ) const
1262 //* Transport the particle on dS, output to P[],C[], for Bz field
1264 const Double_t kCLight = 0.000299792458;
1267 Double_t s = TMath::Sin(bs), c = TMath::Cos(bs);
1269 if( TMath::Abs(bs)>1.e-10){
1273 const Double_t kOvSqr6 = 1./TMath::Sqrt(6.);
1274 sB = (1.-bs*kOvSqr6)*(1.+bs*kOvSqr6)*t;
1278 Double_t px = fP[3];
1279 Double_t py = fP[4];
1280 Double_t pz = fP[5];
1282 p[0] = fP[0] + sB*px + cB*py;
1283 p[1] = fP[1] - cB*px + sB*py;
1284 p[2] = fP[2] + t*pz;
1286 p[4] = -s*px + c*py;
1292 Double_t mJ[8][8] = { {1,0,0, sB, cB, 0, 0, 0 },
1293 {0,1,0, -cB, sB, 0, 0, 0 },
1294 {0,0,1, 0, 0, t, 0, 0 },
1295 {0,0,0, c, s, 0, 0, 0 },
1296 {0,0,0, -s, c, 0, 0, 0 },
1297 {0,0,0, 0, 0, 1, 0, 0 },
1298 {0,0,0, 0, 0, 0, 1, 0 },
1299 {0,0,0, 0, 0, 0, 0, 1 } };
1301 for( Int_t k=0,i=0; i<8; i++)
1302 for( Int_t j=0; j<=i; j++, k++ ) mA[i][j] = mA[j][i] = fC[k];
1305 for( Int_t i=0; i<8; i++ )
1306 for( Int_t j=0; j<8; j++ ){
1308 for( Int_t k=0; k<8; k++ ) mJC[i][j]+=mJ[i][k]*mA[k][j];
1311 for( Int_t k=0,i=0; i<8; i++)
1312 for( Int_t j=0; j<=i; j++, k++ ){
1314 for( Int_t l=0; l<8; l++ ) e[k]+=mJC[i][l]*mJ[j][l];
1321 c6=fC[6], c7=fC[7], c8=fC[8], c17=fC[17], c18=fC[18],
1322 c24 = fC[24], c31 = fC[31];
1326 mJC13 = c7 - cB*fC[9] + sB*fC[13],
1327 mJC14 = fC[11] - cBC13 + sB*fC[14],
1329 mJC24 = fC[12] + t*fC[19],
1330 mJC33 = c*fC[9] + s*fC[13],
1331 mJC34 = c*fC[13] + s*fC[14],
1332 mJC43 = -s*fC[9] + c*fC[13],
1333 mJC44 = -s*fC[13] + c*fC[14];
1336 e[0]= fC[0] + 2*(sB*c6 + cB*fC[10]) + (sB*fC[9] + 2*cBC13)*sB + cB*cB*fC[14];
1337 e[1]= fC[1] - cB*c6 + sB*fC[10] + mJC13*sB + mJC14*cB;
1338 e[2]= fC[2] - cB*c7 + sB*fC[11] - mJC13*cB + mJC14*sB;
1339 e[3]= fC[3] + t*fC[15] + mJC23*sB + mJC24*cB;
1340 e[4]= fC[4] + t*fC[16] - mJC23*cB + mJC24*sB;
1342 e[15]= fC[15] + c18*sB + fC[19]*cB;
1343 e[16]= fC[16] - c18*cB + fC[19]*sB;
1344 e[17]= c17 + fC[20]*t;
1345 e[18]= c18*c + fC[19]*s;
1346 e[19]= -c18*s + fC[19]*c;
1348 e[5]= fC[5] + (c17 + e[17] )*t;
1350 e[6]= c*c6 + s*fC[10] + mJC33*sB + mJC34*cB;
1351 e[7]= c*c7 + s*fC[11] - mJC33*cB + mJC34*sB;
1352 e[8]= c*c8 + s*fC[12] + e[18]*t;
1353 e[9]= mJC33*c + mJC34*s;
1354 e[10]= -s*c6 + c*fC[10] + mJC43*sB + mJC44*cB;
1357 e[11]= -s*c7 + c*fC[11] - mJC43*cB + mJC44*sB;
1358 e[12]= -s*c8 + c*fC[12] + e[19]*t;
1359 e[13]= mJC43*c + mJC44*s;
1360 e[14]= -mJC43*s + mJC44*c;
1362 e[21]= fC[21] + fC[25]*cB + c24*sB;
1363 e[22]= fC[22] - c24*cB + fC[25]*sB;
1364 e[23]= fC[23] + fC[26]*t;
1365 e[24]= c*c24 + s*fC[25];
1366 e[25]= c*fC[25] - c24*s;
1369 e[28]= fC[28] + fC[32]*cB + c31*sB;
1370 e[29]= fC[29] - c31*cB + fC[32]*sB;
1371 e[30]= fC[30] + fC[33]*t;
1372 e[31]= c*c31 + s*fC[32];
1373 e[32]= c*fC[32] - s*c31;
1380 Double_t AliKFParticleBase::GetDistanceFromVertex( const AliKFParticleBase &Vtx ) const
1382 //* Calculate distance from vertex [cm]
1384 return GetDistanceFromVertex( Vtx.fP );
1387 Double_t AliKFParticleBase::GetDistanceFromVertex( const Double_t vtx[] ) const
1389 //* Calculate distance from vertex [cm]
1391 Double_t mP[8], mC[36];
1392 Transport( GetDStoPoint(vtx), mP, mC );
1393 Double_t d[3]={ vtx[0]-mP[0], vtx[1]-mP[1], vtx[2]-mP[2]};
1394 return TMath::Sqrt( d[0]*d[0]+d[1]*d[1]+d[2]*d[2] );
1397 Double_t AliKFParticleBase::GetDistanceFromParticle( const AliKFParticleBase &p )
1400 //* Calculate distance to other particle [cm]
1403 GetDStoParticle( p, dS, dS1 );
1404 Double_t mP[8], mC[36], mP1[8], mC1[36];
1405 Transport( dS, mP, mC );
1406 p.Transport( dS1, mP1, mC1 );
1407 Double_t dx = mP[0]-mP1[0];
1408 Double_t dy = mP[1]-mP1[1];
1409 Double_t dz = mP[2]-mP1[2];
1410 return TMath::Sqrt(dx*dx+dy*dy+dz*dz);
1413 Double_t AliKFParticleBase::GetDeviationFromVertex( const AliKFParticleBase &Vtx ) const
1415 //* Calculate sqrt(Chi2/ndf) deviation from vertex
1417 return GetDeviationFromVertex( Vtx.fP, Vtx.fC );
1421 Double_t AliKFParticleBase::GetDeviationFromVertex( const Double_t v[], const Double_t Cv[] ) const
1423 //* Calculate sqrt(Chi2/ndf) deviation from vertex
1424 //* v = [xyz], Cv=[Cxx,Cxy,Cyy,Cxz,Cyz,Czz]-covariance matrix
1429 Transport( GetDStoPoint(v), mP, mC );
1431 Double_t d[3]={ v[0]-mP[0], v[1]-mP[1], v[2]-mP[2]};
1433 Double_t sigmaS = .1+10.*TMath::Sqrt( (d[0]*d[0]+d[1]*d[1]+d[2]*d[2])/
1434 (mP[3]*mP[3]+mP[4]*mP[4]+mP[5]*mP[5]) );
1437 Double_t h[3] = { mP[3]*sigmaS, mP[4]*sigmaS, mP[5]*sigmaS };
1441 mC[1] +h[1]*h[0], mC[2] +h[1]*h[1],
1442 mC[3] +h[2]*h[0], mC[4] +h[2]*h[1], mC[5] +h[2]*h[2] };
1455 mS[0] = mSi[2]*mSi[5] - mSi[4]*mSi[4];
1456 mS[1] = mSi[3]*mSi[4] - mSi[1]*mSi[5];
1457 mS[2] = mSi[0]*mSi[5] - mSi[3]*mSi[3];
1458 mS[3] = mSi[1]*mSi[4] - mSi[2]*mSi[3];
1459 mS[4] = mSi[1]*mSi[3] - mSi[0]*mSi[4];
1460 mS[5] = mSi[0]*mSi[2] - mSi[1]*mSi[1];
1462 Double_t s = ( mSi[0]*mS[0] + mSi[1]*mS[1] + mSi[3]*mS[3] );
1463 s = ( s > 1.E-20 ) ?1./s :0;
1465 return TMath::Sqrt( TMath::Abs(s*( ( mS[0]*d[0] + mS[1]*d[1] + mS[3]*d[2])*d[0]
1466 +(mS[1]*d[0] + mS[2]*d[1] + mS[4]*d[2])*d[1]
1467 +(mS[3]*d[0] + mS[4]*d[1] + mS[5]*d[2])*d[2] ))/2);
1471 Double_t AliKFParticleBase::GetDeviationFromParticle( const AliKFParticleBase &p )
1474 //* Calculate sqrt(Chi2/ndf) deviation from other particle
1477 GetDStoParticle( p, dS, dS1 );
1478 Double_t mP1[8], mC1[36];
1479 p.Transport( dS1, mP1, mC1 );
1481 Double_t d[3]={ fP[0]-mP1[0], fP[1]-mP1[1], fP[2]-mP1[2]};
1483 Double_t sigmaS = .1+10.*TMath::Sqrt( (d[0]*d[0]+d[1]*d[1]+d[2]*d[2])/
1484 (mP1[3]*mP1[3]+mP1[4]*mP1[4]+mP1[5]*mP1[5]) );
1486 Double_t h[3] = { mP1[3]*sigmaS, mP1[4]*sigmaS, mP1[5]*sigmaS };
1495 return GetDeviationFromVertex( mP1, mC1 )*TMath::Sqrt(2./1.);
1500 void AliKFParticleBase::SubtractFromVertex( AliKFParticleBase &Vtx ) const
1502 //* Subtract the particle from the vertex
1506 GetFieldValue( Vtx.fP, fld );
1507 const Double_t kCLight = 0.000299792458;
1508 fld[0]*=kCLight; fld[1]*=kCLight; fld[2]*=kCLight;
1514 if( Vtx.fIsLinearized ){
1515 GetMeasurement( Vtx.fVtxGuess, m, mCm );
1517 GetMeasurement( Vtx.fP, m, mCm );
1533 Double_t mSi[6] = { mV[0]-Vtx.fC[0],
1534 mV[1]-Vtx.fC[1], mV[2]-Vtx.fC[2],
1535 mV[3]-Vtx.fC[3], mV[4]-Vtx.fC[4], mV[5]-Vtx.fC[5] };
1537 mS[0] = mSi[2]*mSi[5] - mSi[4]*mSi[4];
1538 mS[1] = mSi[3]*mSi[4] - mSi[1]*mSi[5];
1539 mS[2] = mSi[0]*mSi[5] - mSi[3]*mSi[3];
1540 mS[3] = mSi[1]*mSi[4] - mSi[2]*mSi[3];
1541 mS[4] = mSi[1]*mSi[3] - mSi[0]*mSi[4];
1542 mS[5] = mSi[0]*mSi[2] - mSi[1]*mSi[1];
1544 Double_t s = ( mSi[0]*mS[0] + mSi[1]*mS[1] + mSi[3]*mS[3] );
1545 s = ( s > 1.E-20 ) ?1./s :0;
1554 //* Residual (measured - estimated)
1556 Double_t zeta[3] = { m[0]-Vtx.fP[0], m[1]-Vtx.fP[1], m[2]-Vtx.fP[2] };
1558 //* mCHt = mCH' - D'
1560 Double_t mCHt0[3], mCHt1[3], mCHt2[3];
1562 mCHt0[0]=Vtx.fC[ 0] ; mCHt1[0]=Vtx.fC[ 1] ; mCHt2[0]=Vtx.fC[ 3] ;
1563 mCHt0[1]=Vtx.fC[ 1] ; mCHt1[1]=Vtx.fC[ 2] ; mCHt2[1]=Vtx.fC[ 4] ;
1564 mCHt0[2]=Vtx.fC[ 3] ; mCHt1[2]=Vtx.fC[ 4] ; mCHt2[2]=Vtx.fC[ 5] ;
1566 //* Kalman gain K = mCH'*S
1568 Double_t k0[3], k1[3], k2[3];
1570 for(Int_t i=0;i<3;++i){
1571 k0[i] = mCHt0[i]*mS[0] + mCHt1[i]*mS[1] + mCHt2[i]*mS[3];
1572 k1[i] = mCHt0[i]*mS[1] + mCHt1[i]*mS[2] + mCHt2[i]*mS[4];
1573 k2[i] = mCHt0[i]*mS[3] + mCHt1[i]*mS[4] + mCHt2[i]*mS[5];
1576 //* New estimation of the vertex position r += K*zeta
1578 Double_t dChi2 = -(mS[0]*zeta[0] + mS[1]*zeta[1] + mS[3]*zeta[2])*zeta[0]
1579 + (mS[1]*zeta[0] + mS[2]*zeta[1] + mS[4]*zeta[2])*zeta[1]
1580 + (mS[3]*zeta[0] + mS[4]*zeta[1] + mS[5]*zeta[2])*zeta[2];
1582 if( Vtx.fChi2 - dChi2 < 0 ) return;
1584 for(Int_t i=0;i<3;++i)
1585 Vtx.fP[i] -= k0[i]*zeta[0] + k1[i]*zeta[1] + k2[i]*zeta[2];
1587 //* New covariance matrix C -= K*(mCH')'
1589 for(Int_t i=0, k=0;i<3;++i){
1590 for(Int_t j=0;j<=i;++j,++k)
1591 Vtx.fC[k] += k0[i]*mCHt0[j] + k1[i]*mCHt1[j] + k2[i]*mCHt2[j];
1602 void AliKFParticleBase::TransportLine( Double_t dS,
1603 Double_t P[], Double_t C[] ) const
1605 //* Transport the particle as a straight line
1607 P[0] = fP[0] + dS*fP[3];
1608 P[1] = fP[1] + dS*fP[4];
1609 P[2] = fP[2] + dS*fP[5];
1616 Double_t c6 = fC[ 6] + dS*fC[ 9];
1617 Double_t c11 = fC[11] + dS*fC[14];
1618 Double_t c17 = fC[17] + dS*fC[20];
1619 Double_t sc13 = dS*fC[13];
1620 Double_t sc18 = dS*fC[18];
1621 Double_t sc19 = dS*fC[19];
1623 C[ 0] = fC[ 0] + dS*( fC[ 6] + c6 );
1624 C[ 2] = fC[ 2] + dS*( fC[11] + c11 );
1625 C[ 5] = fC[ 5] + dS*( fC[17] + c17 );
1627 C[ 7] = fC[ 7] + sc13;
1628 C[ 8] = fC[ 8] + sc18;
1631 C[12] = fC[12] + sc19;
1633 C[ 1] = fC[ 1] + dS*( fC[10] + C[ 7] );
1634 C[ 3] = fC[ 3] + dS*( fC[15] + C[ 8] );
1635 C[ 4] = fC[ 4] + dS*( fC[16] + C[12] );
1638 C[10] = fC[10] + sc13;
1643 C[15] = fC[15] + sc18;
1644 C[16] = fC[16] + sc19;
1650 C[21] = fC[21] + dS*fC[24];
1651 C[22] = fC[22] + dS*fC[25];
1652 C[23] = fC[23] + dS*fC[26];
1658 C[28] = fC[28] + dS*fC[31];
1659 C[29] = fC[29] + dS*fC[32];
1660 C[30] = fC[30] + dS*fC[33];
1670 void AliKFParticleBase::ConstructGammaBz( const AliKFParticleBase &daughter1,
1671 const AliKFParticleBase &daughter2, double Bz )
1675 const AliKFParticleBase *daughters[2] = { &daughter1, &daughter2};
1679 if( !fIsLinearized ){
1683 daughter1.GetDStoParticle(daughter2, ds, ds1);
1684 daughter1.Transport( ds, m, mCd );
1688 daughter2.Transport( ds1, m, mCd );
1689 v0[0] = .5*( v0[0] + m[0] );
1690 v0[1] = .5*( v0[1] + m[1] );
1691 v0[2] = .5*( v0[2] + m[2] );
1693 v0[0] = fVtxGuess[0];
1694 v0[1] = fVtxGuess[1];
1695 v0[2] = fVtxGuess[2];
1698 fAtProductionVertex = 0;
1710 // fit daughters to the vertex guess
1712 double daughterP[2][8], daughterC[2][36];
1713 double vtxMom[2][3];
1716 for( int id=0; id<2; id++ ){
1718 double *p = daughterP[id];
1719 double *mC = daughterC[id];
1721 daughters[id]->GetMeasurement( v0, p, mC );
1724 InvertSym3(mC, mAi );
1728 mB[0][0] = mC[ 6]*mAi[0] + mC[ 7]*mAi[1] + mC[ 8]*mAi[3];
1729 mB[0][1] = mC[ 6]*mAi[1] + mC[ 7]*mAi[2] + mC[ 8]*mAi[4];
1730 mB[0][2] = mC[ 6]*mAi[3] + mC[ 7]*mAi[4] + mC[ 8]*mAi[5];
1732 mB[1][0] = mC[10]*mAi[0] + mC[11]*mAi[1] + mC[12]*mAi[3];
1733 mB[1][1] = mC[10]*mAi[1] + mC[11]*mAi[2] + mC[12]*mAi[4];
1734 mB[1][2] = mC[10]*mAi[3] + mC[11]*mAi[4] + mC[12]*mAi[5];
1736 mB[2][0] = mC[15]*mAi[0] + mC[16]*mAi[1] + mC[17]*mAi[3];
1737 mB[2][1] = mC[15]*mAi[1] + mC[16]*mAi[2] + mC[17]*mAi[4];
1738 mB[2][2] = mC[15]*mAi[3] + mC[16]*mAi[4] + mC[17]*mAi[5];
1740 Double_t z[3] = { v0[0]-p[0], v0[1]-p[1], v0[2]-p[2] };
1742 vtxMom[id][0] = p[3] + mB[0][0]*z[0] + mB[0][1]*z[1] + mB[0][2]*z[2];
1743 vtxMom[id][1] = p[4] + mB[1][0]*z[0] + mB[1][1]*z[1] + mB[1][2]*z[2];
1744 vtxMom[id][2] = p[5] + mB[2][0]*z[0] + mB[2][1]*z[1] + mB[2][2]*z[2];
1746 daughters[id]->Transport( daughters[id]->GetDStoPoint(v0), p, mC );
1750 } // fit daughters to guess
1756 double mpx0 = vtxMom[0][0]+vtxMom[1][0];
1757 double mpy0 = vtxMom[0][1]+vtxMom[1][1];
1758 double mpt0 = TMath::Sqrt(mpx0*mpx0 + mpy0*mpy0);
1759 // double a0 = TMath::ATan2(mpy0,mpx0);
1761 double ca0 = mpx0/mpt0;
1762 double sa0 = mpy0/mpt0;
1763 double r[3] = { v0[0], v0[1], v0[2] };
1764 double mC[3][3] = {{10000., 0 , 0 },
1769 for( int id=0; id<2; id++ ){
1770 const Double_t kCLight = 0.000299792458;
1771 Double_t q = Bz*daughters[id]->GetQ()*kCLight;
1772 Double_t px0 = vtxMom[id][0];
1773 Double_t py0 = vtxMom[id][1];
1774 Double_t pz0 = vtxMom[id][2];
1775 Double_t pt0 = TMath::Sqrt(px0*px0+py0*py0);
1776 Double_t mG[3][6], mB[3], mH[3][3];
1778 // m = {x,y,z,Px,Py,Pz};
1781 // q*x + Py - q*vx - sin(a)*Pt = 0
1782 // q*y - Px - q*vy + cos(a)*Pt = 0
1783 // (Px*cos(a) + Py*sin(a) ) (vz -z) - Pz( cos(a)*(vx-x) + sin(a)*(vy-y)) = 0
1788 mG[0][3] = -sa0*px0/pt0;
1789 mG[0][4] = 1 -sa0*py0/pt0;
1794 mB[0] = py0 - sa0*pt0 - mG[0][3]*px0 - mG[0][4]*py0 ;
1796 // q*y - Px - q*vy + cos(a)*Pt = 0
1801 mG[1][3] = -1 + ca0*px0/pt0;
1802 mG[1][4] = + ca0*py0/pt0;
1807 mB[1] = -px0 + ca0*pt0 - mG[1][3]*px0 - mG[1][4]*py0 ;
1809 // (Px*cos(a) + Py*sin(a) ) (z -vz) - Pz( cos(a)*(x-vx) + sin(a)*(y-vy)) = 0
1811 mG[2][0] = -pz0*ca0;
1812 mG[2][1] = -pz0*sa0;
1813 mG[2][2] = px0*ca0 + py0*sa0;
1818 mH[2][0] = mG[2][0];
1819 mH[2][1] = mG[2][1];
1820 mH[2][2] = mG[2][2];
1831 for( int i=0; i<3; i++ ){
1833 for( int k=0; k<6; k++ ) m[i]+=mG[i][k]*daughterP[id][k];
1835 for( int i=0; i<3; i++ ){
1836 for( int j=0; j<6; j++ ){
1838 for( int k=0; k<6; k++ ) mGV[i][j]+=mG[i][k]*daughterC[id][ IJ(k,j) ];
1841 for( int i=0, k=0; i<3; i++ ){
1842 for( int j=0; j<=i; j++,k++ ){
1844 for( int l=0; l<6; l++ ) mV[k]+=mGV[i][l]*mG[j][l];
1851 Double_t mCHt[3][3];
1854 for( int i=0; i<3; i++ ){
1856 for( int k=0; k<3; k++ ) mHr[i]+= mH[i][k]*r[k];
1859 for( int i=0; i<3; i++ ){
1860 for( int j=0; j<3; j++){
1862 for( int k=0; k<3; k++ ) mCHt[i][j]+= mC[i][k]*mH[j][k];
1866 for( int i=0, k=0; i<3; i++ ){
1867 for( int j=0; j<=i; j++, k++ ){
1869 for( int l=0; l<3; l++ ) mHCHt[k]+= mH[i][l]*mCHt[l][j];
1873 Double_t mS[6] = { mHCHt[0]+mV[0],
1874 mHCHt[1]+mV[1], mHCHt[2]+mV[2],
1875 mHCHt[3]+mV[3], mHCHt[4]+mV[4], mHCHt[5]+mV[5] };
1880 //* Residual (measured - estimated)
1882 Double_t zeta[3] = { m[0]-mHr[0], m[1]-mHr[1], m[2]-mHr[2] };
1884 //* Kalman gain K = mCH'*S
1888 for(Int_t i=0;i<3;++i){
1889 k[i][0] = mCHt[i][0]*mS[0] + mCHt[i][1]*mS[1] + mCHt[i][2]*mS[3];
1890 k[i][1] = mCHt[i][0]*mS[1] + mCHt[i][1]*mS[2] + mCHt[i][2]*mS[4];
1891 k[i][2] = mCHt[i][0]*mS[3] + mCHt[i][1]*mS[4] + mCHt[i][2]*mS[5];
1894 //* New estimation of the vertex position r += K*zeta
1896 for(Int_t i=0;i<3;++i)
1897 r[i] = r[i] + k[i][0]*zeta[0] + k[i][1]*zeta[1] + k[i][2]*zeta[2];
1899 //* New covariance matrix C -= K*(mCH')'
1901 for(Int_t i=0;i<3;++i){
1902 for(Int_t j=0;j<=i;++j){
1903 mC[i][j] = mC[i][j] - (k[i][0]*mCHt[j][0] + k[i][1]*mCHt[j][1] + k[i][2]*mCHt[j][2]);
1904 mC[j][i] = mC[i][j];
1911 chi2 += ( ( mS[0]*zeta[0] + mS[1]*zeta[1] + mS[3]*zeta[2] )*zeta[0]
1912 +(mS[1]*zeta[0] + mS[2]*zeta[1] + mS[4]*zeta[2] )*zeta[1]
1913 +(mS[3]*zeta[0] + mS[4]*zeta[1] + mS[5]*zeta[2] )*zeta[2] );
1920 for( int i=0; i<3; i++ ) fP[i] = r[i];
1921 for( int i=0,k=0; i<3; i++ ){
1922 for( int j=0; j<=i; j++,k++ ){
1928 // now fit daughters to the vertex
1933 for(Int_t i=3;i<8;++i) fP[i]=0.;
1934 for(Int_t i=6;i<36;++i) fC[i]=0.;
1937 for( int id=0; id<2; id++ ){
1939 double *p = daughterP[id];
1940 double *mC = daughterC[id];
1941 daughters[id]->GetMeasurement( v0, p, mC );
1943 const Double_t *m = fP, *mV = fC;
1946 InvertSym3(mC, mAi );
1950 mB[0][0] = mC[ 6]*mAi[0] + mC[ 7]*mAi[1] + mC[ 8]*mAi[3];
1951 mB[0][1] = mC[ 6]*mAi[1] + mC[ 7]*mAi[2] + mC[ 8]*mAi[4];
1952 mB[0][2] = mC[ 6]*mAi[3] + mC[ 7]*mAi[4] + mC[ 8]*mAi[5];
1954 mB[1][0] = mC[10]*mAi[0] + mC[11]*mAi[1] + mC[12]*mAi[3];
1955 mB[1][1] = mC[10]*mAi[1] + mC[11]*mAi[2] + mC[12]*mAi[4];
1956 mB[1][2] = mC[10]*mAi[3] + mC[11]*mAi[4] + mC[12]*mAi[5];
1958 mB[2][0] = mC[15]*mAi[0] + mC[16]*mAi[1] + mC[17]*mAi[3];
1959 mB[2][1] = mC[15]*mAi[1] + mC[16]*mAi[2] + mC[17]*mAi[4];
1960 mB[2][2] = mC[15]*mAi[3] + mC[16]*mAi[4] + mC[17]*mAi[5];
1962 mB[3][0] = mC[21]*mAi[0] + mC[22]*mAi[1] + mC[23]*mAi[3];
1963 mB[3][1] = mC[21]*mAi[1] + mC[22]*mAi[2] + mC[23]*mAi[4];
1964 mB[3][2] = mC[21]*mAi[3] + mC[22]*mAi[4] + mC[23]*mAi[5];
1967 Double_t z[3] = { m[0]-p[0], m[1]-p[1], m[2]-p[2] };
1970 Double_t mAV[6] = { mC[0]-mV[0], mC[1]-mV[1], mC[2]-mV[2],
1971 mC[3]-mV[3], mC[4]-mV[4], mC[5]-mV[5] };
1974 if( !InvertSym3(mAV, mAVi) ){
1975 Double_t dChi2 = ( +(mAVi[0]*z[0] + mAVi[1]*z[1] + mAVi[3]*z[2])*z[0]
1976 +(mAVi[1]*z[0] + mAVi[2]*z[1] + mAVi[4]*z[2])*z[1]
1977 +(mAVi[3]*z[0] + mAVi[4]*z[1] + mAVi[5]*z[2])*z[2] );
1978 fChi2+= TMath::Abs( dChi2 );
1983 //* Add the daughter momentum to the particle momentum
1985 fP[3]+= p[3] + mB[0][0]*z[0] + mB[0][1]*z[1] + mB[0][2]*z[2];
1986 fP[4]+= p[4] + mB[1][0]*z[0] + mB[1][1]*z[1] + mB[1][2]*z[2];
1987 fP[5]+= p[5] + mB[2][0]*z[0] + mB[2][1]*z[1] + mB[2][2]*z[2];
1988 fP[6]+= p[6] + mB[3][0]*z[0] + mB[3][1]*z[1] + mB[3][2]*z[2];
1990 Double_t d0, d1, d2;
1992 d0= mB[0][0]*mV[0] + mB[0][1]*mV[1] + mB[0][2]*mV[3] - mC[ 6];
1993 d1= mB[0][0]*mV[1] + mB[0][1]*mV[2] + mB[0][2]*mV[4] - mC[ 7];
1994 d2= mB[0][0]*mV[3] + mB[0][1]*mV[4] + mB[0][2]*mV[5] - mC[ 8];
1996 //fC[6]+= mC[ 6] + d0;
1997 //fC[7]+= mC[ 7] + d1;
1998 //fC[8]+= mC[ 8] + d2;
1999 fC[9]+= mC[ 9] + d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
2001 d0= mB[1][0]*mV[0] + mB[1][1]*mV[1] + mB[1][2]*mV[3] - mC[10];
2002 d1= mB[1][0]*mV[1] + mB[1][1]*mV[2] + mB[1][2]*mV[4] - mC[11];
2003 d2= mB[1][0]*mV[3] + mB[1][1]*mV[4] + mB[1][2]*mV[5] - mC[12];
2005 //fC[10]+= mC[10]+ d0;
2006 //fC[11]+= mC[11]+ d1;
2007 //fC[12]+= mC[12]+ d2;
2008 fC[13]+= mC[13]+ d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
2009 fC[14]+= mC[14]+ d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
2011 d0= mB[2][0]*mV[0] + mB[2][1]*mV[1] + mB[2][2]*mV[3] - mC[15];
2012 d1= mB[2][0]*mV[1] + mB[2][1]*mV[2] + mB[2][2]*mV[4] - mC[16];
2013 d2= mB[2][0]*mV[3] + mB[2][1]*mV[4] + mB[2][2]*mV[5] - mC[17];
2015 //fC[15]+= mC[15]+ d0;
2016 //fC[16]+= mC[16]+ d1;
2017 //fC[17]+= mC[17]+ d2;
2018 fC[18]+= mC[18]+ d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
2019 fC[19]+= mC[19]+ d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
2020 fC[20]+= mC[20]+ d0*mB[2][0] + d1*mB[2][1] + d2*mB[2][2];
2022 d0= mB[3][0]*mV[0] + mB[3][1]*mV[1] + mB[3][2]*mV[3] - mC[21];
2023 d1= mB[3][0]*mV[1] + mB[3][1]*mV[2] + mB[3][2]*mV[4] - mC[22];
2024 d2= mB[3][0]*mV[3] + mB[3][1]*mV[4] + mB[3][2]*mV[5] - mC[23];
2026 //fC[21]+= mC[21] + d0;
2027 //fC[22]+= mC[22] + d1;
2028 //fC[23]+= mC[23] + d2;
2029 fC[24]+= mC[24] + d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
2030 fC[25]+= mC[25] + d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
2031 fC[26]+= mC[26] + d0*mB[2][0] + d1*mB[2][1] + d2*mB[2][2];
2032 fC[27]+= mC[27] + d0*mB[3][0] + d1*mB[3][1] + d2*mB[3][2];
2035 SetMassConstraint(0,0);
2039 Bool_t AliKFParticleBase::InvertSym3( const Double_t A[], Double_t Ai[] )
2041 //* Invert symmetric matric stored in low-triagonal form
2044 double a0 = A[0], a1 = A[1], a2 = A[2], a3 = A[3];
2046 Ai[0] = a2*A[5] - A[4]*A[4];
2047 Ai[1] = a3*A[4] - a1*A[5];
2048 Ai[3] = a1*A[4] - a2*a3;
2049 Double_t det = (a0*Ai[0] + a1*Ai[1] + a3*Ai[3]);
2050 if( det>1.e-20 ) det = 1./det;
2058 Ai[2] = ( a0*A[5] - a3*a3 )*det;
2059 Ai[4] = ( a1*a3 - a0*A[4] )*det;
2060 Ai[5] = ( a0*a2 - a1*a1 )*det;
2064 void AliKFParticleBase::MultQSQt( const Double_t Q[], const Double_t S[], Double_t SOut[] )
2066 //* Matrix multiplication Q*S*Q^T, Q - square matrix, S - symmetric
2071 for( Int_t i=0, ij=0; i<kN; i++ ){
2072 for( Int_t j=0; j<kN; j++, ++ij ){
2074 for( Int_t k=0; k<kN; ++k ) mA[ij]+= S[( k<=i ) ? i*(i+1)/2+k :k*(k+1)/2+i] * Q[ j*kN+k];
2078 for( Int_t i=0; i<kN; i++ ){
2079 for( Int_t j=0; j<=i; j++ ){
2080 Int_t ij = ( j<=i ) ? i*(i+1)/2+j :j*(j+1)/2+i;
2082 for( Int_t k=0; k<kN; k++ ) SOut[ij] += Q[ i*kN+k ] * mA[ k*kN+j ];
2088 // 72-charachters line to define the printer border
2089 //3456789012345678901234567890123456789012345678901234567890123456789012