ROOT logo
//----------------------------------------------------------------------------
// Implementation of the AliKFParticle class
// .
// @author  S.Gorbunov, I.Kisel
// @version 1.0
// @since   13.05.07
// 
// Class to reconstruct and store the decayed particle parameters.
// The method is described in CBM-SOFT note 2007-003, 
// ``Reconstruction of decayed particles based on the Kalman filter'', 
// http://www.gsi.de/documents/DOC-2007-May-14-1.pdf
//
// This class is ALICE interface to general mathematics in AliKFParticleCore
// 
//  -= Copyright &copy ALICE HLT Group =-
//____________________________________________________________________________


#include "AliKFParticle.h"
#include "TDatabasePDG.h"
#include "TParticlePDG.h"
#include "AliVTrack.h"
#include "AliVVertex.h"

#include "AliExternalTrackParam.h"

ClassImp(AliKFParticle)

Double_t AliKFParticle::fgBz = -5.;  //* Bz compoment of the magnetic field

AliKFParticle::AliKFParticle( const AliKFParticle &d1, const AliKFParticle &d2, Bool_t gamma )
{
  if (!gamma) {
    AliKFParticle mother;
    mother+= d1;
    mother+= d2;
    *this = mother;
  } else
    ConstructGamma(d1, d2);
}

void AliKFParticle::Create( const Double_t Param[], const Double_t Cov[], Int_t Charge, Int_t PID )
{
  // Constructor from "cartesian" track, PID hypothesis should be provided
  //
  // Param[6] = { X, Y, Z, Px, Py, Pz } - position and momentum
  // Cov [21] = lower-triangular part of the covariance matrix:
  //
  //                (  0  .  .  .  .  . )
  //                (  1  2  .  .  .  . )
  //  Cov. matrix = (  3  4  5  .  .  . ) - numbering of covariance elements in Cov[]
  //                (  6  7  8  9  .  . )
  //                ( 10 11 12 13 14  . )
  //                ( 15 16 17 18 19 20 )
  Double_t C[21];
  for( int i=0; i<21; i++ ) C[i] = Cov[i];
  
  TParticlePDG* particlePDG = TDatabasePDG::Instance()->GetParticle(PID);
  Double_t mass = (particlePDG) ? particlePDG->Mass() :0.13957;
  
  AliKFParticleBase::Initialize( Param, C, Charge, mass );
}

void AliKFParticle::Create( const Double_t Param[], const Double_t Cov[], Int_t Charge, Double_t Mass )
{
  // Constructor from "cartesian" track, PID hypothesis should be provided
  //
  // Param[6] = { X, Y, Z, Px, Py, Pz } - position and momentum
  // Cov [21] = lower-triangular part of the covariance matrix:
  //
  //                (  0  .  .  .  .  . )
  //                (  1  2  .  .  .  . )
  //  Cov. matrix = (  3  4  5  .  .  . ) - numbering of covariance elements in Cov[]
  //                (  6  7  8  9  .  . )
  //                ( 10 11 12 13 14  . )
  //                ( 15 16 17 18 19 20 )
  Double_t C[21];
  for( int i=0; i<21; i++ ) C[i] = Cov[i];

  AliKFParticleBase::Initialize( Param, C, Charge, Mass );
}

AliKFParticle::AliKFParticle( const AliVTrack &track, Int_t PID )
{
  // Constructor from ALICE track, PID hypothesis should be provided

  track.GetXYZ(fP);
  track.PxPyPz(fP+3);
  fQ = track.Charge();
  track.GetCovarianceXYZPxPyPz( fC );
  Create(fP,fC,fQ,PID);
}

AliKFParticle::AliKFParticle( const AliExternalTrackParam &track, Double_t Mass, Int_t Charge )
{
  // Constructor from ALICE track, Mass and Charge hypothesis should be provided

  track.GetXYZ(fP);
  track.PxPyPz(fP+3);
  fQ = track.Charge()*TMath::Abs(Charge);
  fP[3] *= TMath::Abs(Charge);
  fP[4] *= TMath::Abs(Charge);
  fP[5] *= TMath::Abs(Charge);

  Double_t pt=1./TMath::Abs(track.GetParameter()[4]) * TMath::Abs(Charge);
  Double_t cs=TMath::Cos(track.GetAlpha()), sn=TMath::Sin(track.GetAlpha());
  Double_t r=TMath::Sqrt((1.-track.GetParameter()[2])*(1.+track.GetParameter()[2]));

  Double_t m00=-sn, m10=cs;
  Double_t m23=-pt*(sn + track.GetParameter()[2]*cs/r), m43=-pt*pt*(r*cs - track.GetParameter()[2]*sn);
  Double_t m24= pt*(cs - track.GetParameter()[2]*sn/r), m44=-pt*pt*(r*sn + track.GetParameter()[2]*cs);
  Double_t m35=pt, m45=-pt*pt*track.GetParameter()[3];

  m43*=track.GetSign();
  m44*=track.GetSign();
  m45*=track.GetSign();

  const Double_t *cTr = track.GetCovariance();

  fC[0 ] = cTr[0]*m00*m00;
  fC[1 ] = cTr[0]*m00*m10; 
  fC[2 ] = cTr[0]*m10*m10;
  fC[3 ] = cTr[1]*m00; 
  fC[4 ] = cTr[1]*m10; 
  fC[5 ] = cTr[2];
  fC[6 ] = m00*(cTr[3]*m23 + cTr[10]*m43); 
  fC[7 ] = m10*(cTr[3]*m23 + cTr[10]*m43); 
  fC[8 ] = cTr[4]*m23 + cTr[11]*m43; 
  fC[9 ] = m23*(cTr[5]*m23 + cTr[12]*m43)  +  m43*(cTr[12]*m23 + cTr[14]*m43);
  fC[10] = m00*(cTr[3]*m24 + cTr[10]*m44); 
  fC[11] = m10*(cTr[3]*m24 + cTr[10]*m44); 
  fC[12] = cTr[4]*m24 + cTr[11]*m44; 
  fC[13] = m23*(cTr[5]*m24 + cTr[12]*m44)  +  m43*(cTr[12]*m24 + cTr[14]*m44);
  fC[14] = m24*(cTr[5]*m24 + cTr[12]*m44)  +  m44*(cTr[12]*m24 + cTr[14]*m44);
  fC[15] = m00*(cTr[6]*m35 + cTr[10]*m45); 
  fC[16] = m10*(cTr[6]*m35 + cTr[10]*m45); 
  fC[17] = cTr[7]*m35 + cTr[11]*m45; 
  fC[18] = m23*(cTr[8]*m35 + cTr[12]*m45)  +  m43*(cTr[13]*m35 + cTr[14]*m45);
  fC[19] = m24*(cTr[8]*m35 + cTr[12]*m45)  +  m44*(cTr[13]*m35 + cTr[14]*m45); 
  fC[20] = m35*(cTr[9]*m35 + cTr[13]*m45)  +  m45*(cTr[13]*m35 + cTr[14]*m45);

  Create(fP,fC,fQ,Mass);
}

AliKFParticle::AliKFParticle( const AliVVertex &vertex )
{
  // Constructor from ALICE vertex

  vertex.GetXYZ( fP );
  vertex.GetCovarianceMatrix( fC );  
  fChi2 = vertex.GetChi2();
  fNDF = 2*vertex.GetNContributors() - 3;
  fQ = 0;
  fAtProductionVertex = 0;
  fIsLinearized = 0;
  fSFromDecay = 0;
}

void AliKFParticle::GetExternalTrackParam( const AliKFParticleBase &p, Double_t &X, Double_t &Alpha, Double_t P[5] ) 
{
  // Conversion to AliExternalTrackParam parameterization

  Double_t cosA = p.GetPx(), sinA = p.GetPy(); 
  Double_t pt = TMath::Sqrt(cosA*cosA + sinA*sinA);
  Double_t pti = 0;
  if( pt<1.e-4 ){
    cosA = 1;
    sinA = 0;
  } else {
    pti = 1./pt;
    cosA*=pti;
    sinA*=pti;
  }
  Alpha = TMath::ATan2(sinA,cosA);  
  X   = p.GetX()*cosA + p.GetY()*sinA;   
  P[0]= p.GetY()*cosA - p.GetX()*sinA;
  P[1]= p.GetZ();
  P[2]= 0;
  P[3]= p.GetPz()*pti;
  P[4]= p.GetQ()*pti;
}

Bool_t AliKFParticle::GetDistanceFromVertexXY( const Double_t vtx[], const Double_t Cv[], Double_t &val, Double_t &err ) const
{
  //* Calculate DCA distance from vertex (transverse impact parameter) in XY
  //* v = [xy], Cv=[Cxx,Cxy,Cyy ]-covariance matrix

  Bool_t ret = 0;
  
  Double_t mP[8];
  Double_t mC[36];
  
  Transport( GetDStoPoint(vtx), mP, mC );  

  Double_t dx = mP[0] - vtx[0];
  Double_t dy = mP[1] - vtx[1];
  Double_t px = mP[3];
  Double_t py = mP[4];
  Double_t pt = TMath::Sqrt(px*px + py*py);
  Double_t ex=0, ey=0;
  if( pt<1.e-4 ){
    ret = 1;
    pt = 1.;
    val = 1.e4;
  } else{
    ex = px/pt;
    ey = py/pt;
    val = dy*ex - dx*ey;
  }

  Double_t h0 = -ey;
  Double_t h1 = ex;
  Double_t h3 = (dy*ey + dx*ex)*ey/pt;
  Double_t h4 = -(dy*ey + dx*ex)*ex/pt;
  
  err = 
    h0*(h0*GetCovariance(0,0) + h1*GetCovariance(0,1) + h3*GetCovariance(0,3) + h4*GetCovariance(0,4) ) +
    h1*(h0*GetCovariance(1,0) + h1*GetCovariance(1,1) + h3*GetCovariance(1,3) + h4*GetCovariance(1,4) ) +
    h3*(h0*GetCovariance(3,0) + h1*GetCovariance(3,1) + h3*GetCovariance(3,3) + h4*GetCovariance(3,4) ) +
    h4*(h0*GetCovariance(4,0) + h1*GetCovariance(4,1) + h3*GetCovariance(4,3) + h4*GetCovariance(4,4) );

  if( Cv ){
    err+= h0*(h0*Cv[0] + h1*Cv[1] ) + h1*(h0*Cv[1] + h1*Cv[2] ); 
  }

  err = TMath::Sqrt(TMath::Abs(err));

  return ret;
}

Bool_t AliKFParticle::GetDistanceFromVertexXY( const Double_t vtx[], Double_t &val, Double_t &err ) const
{
  return GetDistanceFromVertexXY( vtx, 0, val, err );
}


Bool_t AliKFParticle::GetDistanceFromVertexXY( const AliKFParticle &Vtx, Double_t &val, Double_t &err ) const 
{
  //* Calculate distance from vertex [cm] in XY-plane

  return GetDistanceFromVertexXY( Vtx.fP, Vtx.fC, val, err );
}

Bool_t AliKFParticle::GetDistanceFromVertexXY( const AliVVertex &Vtx, Double_t &val, Double_t &err ) const 
{
  //* Calculate distance from vertex [cm] in XY-plane

  return GetDistanceFromVertexXY( AliKFParticle(Vtx), val, err );
}

Double_t AliKFParticle::GetDistanceFromVertexXY( const Double_t vtx[] ) const
{
  //* Calculate distance from vertex [cm] in XY-plane
  Double_t val, err;
  GetDistanceFromVertexXY( vtx, 0, val, err );
  return val;
}

Double_t AliKFParticle::GetDistanceFromVertexXY( const AliKFParticle &Vtx ) const 
{
  //* Calculate distance from vertex [cm] in XY-plane

  return GetDistanceFromVertexXY( Vtx.fP );
}

Double_t AliKFParticle::GetDistanceFromVertexXY( const AliVVertex &Vtx ) const 
{
  //* Calculate distance from vertex [cm] in XY-plane

  return GetDistanceFromVertexXY( AliKFParticle(Vtx).fP );
}

Double_t AliKFParticle::GetDistanceFromParticleXY( const AliKFParticle &p ) const 
{
  //* Calculate distance to other particle [cm]

  Double_t dS, dS1;
  GetDStoParticleXY( p, dS, dS1 );   
  Double_t mP[8], mC[36], mP1[8], mC1[36];
  Transport( dS, mP, mC ); 
  p.Transport( dS1, mP1, mC1 ); 
  Double_t dx = mP[0]-mP1[0]; 
  Double_t dy = mP[1]-mP1[1]; 
  return TMath::Sqrt(dx*dx+dy*dy);
}

Double_t AliKFParticle::GetDeviationFromParticleXY( const AliKFParticle &p ) const 
{
  //* Calculate sqrt(Chi2/ndf) deviation from other particle

  Double_t dS, dS1;
  GetDStoParticleXY( p, dS, dS1 );   
  Double_t mP1[8], mC1[36];
  p.Transport( dS1, mP1, mC1 ); 

  Double_t d[2]={ fP[0]-mP1[0], fP[1]-mP1[1] };

  Double_t sigmaS = .1+10.*TMath::Sqrt( (d[0]*d[0]+d[1]*d[1] )/
					(mP1[3]*mP1[3]+mP1[4]*mP1[4] )  );

  Double_t h[2] = { mP1[3]*sigmaS, mP1[4]*sigmaS };       
  
  mC1[0] +=h[0]*h[0];
  mC1[1] +=h[1]*h[0]; 
  mC1[2] +=h[1]*h[1]; 

  return GetDeviationFromVertexXY( mP1, mC1 )*TMath::Sqrt(2./1.);
}


Double_t AliKFParticle::GetDeviationFromVertexXY( const Double_t vtx[], const Double_t Cv[] ) const 
{
  //* Calculate sqrt(Chi2/ndf) deviation from vertex
  //* v = [xyz], Cv=[Cxx,Cxy,Cyy,Cxz,Cyz,Czz]-covariance matrix

  Double_t val, err;
  Bool_t problem = GetDistanceFromVertexXY( vtx, Cv, val, err );
  if( problem || err<1.e-20 ) return 1.e4;
  else return val/err;
}


Double_t AliKFParticle::GetDeviationFromVertexXY( const AliKFParticle &Vtx ) const  
{
  //* Calculate sqrt(Chi2/ndf) deviation from vertex
  //* v = [xyz], Cv=[Cxx,Cxy,Cyy,Cxz,Cyz,Czz]-covariance matrix

  return GetDeviationFromVertexXY( Vtx.fP, Vtx.fC );
}

Double_t AliKFParticle::GetDeviationFromVertexXY( const AliVVertex &Vtx ) const 
{
  //* Calculate sqrt(Chi2/ndf) deviation from vertex
  //* v = [xyz], Cv=[Cxx,Cxy,Cyy,Cxz,Cyz,Czz]-covariance matrix

  AliKFParticle v(Vtx);
  return GetDeviationFromVertexXY( v.fP, v.fC );
}

Double_t AliKFParticle::GetAngle  ( const AliKFParticle &p ) const 
{
  //* Calculate the opening angle between two particles

  Double_t dS, dS1;
  GetDStoParticle( p, dS, dS1 );   
  Double_t mP[8], mC[36], mP1[8], mC1[36];
  Transport( dS, mP, mC ); 
  p.Transport( dS1, mP1, mC1 ); 
  Double_t n = TMath::Sqrt( mP[3]*mP[3] + mP[4]*mP[4] + mP[5]*mP[5] );
  Double_t n1= TMath::Sqrt( mP1[3]*mP1[3] + mP1[4]*mP1[4] + mP1[5]*mP1[5] );
  n*=n1;
  Double_t a = 0;
  if( n>1.e-8 ) a = ( mP[3]*mP1[3] + mP[4]*mP1[4] + mP[5]*mP1[5] )/n;
  if (TMath::Abs(a)<1.) a = TMath::ACos(a);
  else a = (a>=0) ?0 :TMath::Pi();
  return a;
}

Double_t AliKFParticle::GetAngleXY( const AliKFParticle &p ) const 
{
  //* Calculate the opening angle between two particles in XY plane

  Double_t dS, dS1;
  GetDStoParticleXY( p, dS, dS1 );   
  Double_t mP[8], mC[36], mP1[8], mC1[36];
  Transport( dS, mP, mC ); 
  p.Transport( dS1, mP1, mC1 ); 
  Double_t n = TMath::Sqrt( mP[3]*mP[3] + mP[4]*mP[4] );
  Double_t n1= TMath::Sqrt( mP1[3]*mP1[3] + mP1[4]*mP1[4] );
  n*=n1;
  Double_t a = 0;
  if( n>1.e-8 ) a = ( mP[3]*mP1[3] + mP[4]*mP1[4] )/n;
  if (TMath::Abs(a)<1.) a = TMath::ACos(a);
  else a = (a>=0) ?0 :TMath::Pi();
  return a;
}

Double_t AliKFParticle::GetAngleRZ( const AliKFParticle &p ) const 
{
  //* Calculate the opening angle between two particles in RZ plane

  Double_t dS, dS1;
  GetDStoParticle( p, dS, dS1 );   
  Double_t mP[8], mC[36], mP1[8], mC1[36];
  Transport( dS, mP, mC ); 
  p.Transport( dS1, mP1, mC1 ); 
  Double_t nr = TMath::Sqrt( mP[3]*mP[3] + mP[4]*mP[4] );
  Double_t n1r= TMath::Sqrt( mP1[3]*mP1[3] + mP1[4]*mP1[4]  );
  Double_t n = TMath::Sqrt( nr*nr + mP[5]*mP[5] );
  Double_t n1= TMath::Sqrt( n1r*n1r + mP1[5]*mP1[5] );
  n*=n1;
  Double_t a = 0;
  if( n>1.e-8 ) a = ( nr*n1r +mP[5]*mP1[5])/n; 
  if (TMath::Abs(a)<1.) a = TMath::ACos(a);
  else a = (a>=0) ?0 :TMath::Pi();
  return a;
}


/*

#include "AliExternalTrackParam.h"

void AliKFParticle::GetDStoParticleALICE( const AliKFParticleBase &p, 
					   Double_t &DS, Double_t &DS1 ) 
  const
{ 
  DS = DS1 = 0;   
  Double_t x1, a1, x2, a2;
  Double_t par1[5], par2[5], cov[15];
  for(int i=0; i<15; i++) cov[i] = 0;
  cov[0] = cov[2] = cov[5] = cov[9] = cov[14] = .001;

  GetExternalTrackParam( *this, x1, a1, par1 );
  GetExternalTrackParam( p, x2, a2, par2 );

  AliExternalTrackParam t1(x1,a1, par1, cov);
  AliExternalTrackParam t2(x2,a2, par2, cov);

  Double_t xe1=0, xe2=0;
  t1.GetDCA( &t2, -GetFieldAlice(), xe1, xe2 );
  t1.PropagateTo( xe1, -GetFieldAlice() );
  t2.PropagateTo( xe2, -GetFieldAlice() );

  Double_t xyz1[3], xyz2[3];
  t1.GetXYZ( xyz1 );
  t2.GetXYZ( xyz2 );
  
  DS = GetDStoPoint( xyz1 );
  DS1 = p.GetDStoPoint( xyz2 );

  return;
}
*/

  // * Pseudo Proper Time of decay = (r*pt) / |pt| * M/|pt|
Double_t AliKFParticle::GetPseudoProperDecayTime( const AliKFParticle &pV, const Double_t& mass, Double_t* timeErr2 ) const
{ // TODO optimize with respect to time and stability
  const Double_t ipt2 = 1/( Px()*Px() + Py()*Py() );
  const Double_t mipt2 = mass*ipt2;
  const Double_t dx = X() - pV.X();
  const Double_t dy = Y() - pV.Y();

  if ( timeErr2 ) {
      // -- calculate error = sigma(f(r)) = f'Cf'
      // r = {x,y,px,py,x_pV,y_pV}
      // df/dr = { px*m/pt^2,
      //     py*m/pt^2,
      //    ( x - x_pV )*m*(1/pt^2 - 2(px/pt^2)^2),
      //    ( y - y_pV )*m*(1/pt^2 - 2(py/pt^2)^2),
      //     -px*m/pt^2,
      //     -py*m/pt^2 }
    const Double_t f0 = Px()*mipt2;
    const Double_t f1 = Py()*mipt2;
    const Double_t mipt2derivative = mipt2*(1-2*Px()*Px()*ipt2);
    const Double_t f2 = dx*mipt2derivative;
    const Double_t f3 = -dy*mipt2derivative;
    const Double_t f4 = -f0;
    const Double_t f5 = -f1;

    const Double_t& mC00 =    GetCovariance(0,0);
    const Double_t& mC10 =    GetCovariance(0,1);
    const Double_t& mC11 =    GetCovariance(1,1);
    const Double_t& mC20 =    GetCovariance(3,0);
    const Double_t& mC21 =    GetCovariance(3,1);
    const Double_t& mC22 =    GetCovariance(3,3);
    const Double_t& mC30 =    GetCovariance(4,0);
    const Double_t& mC31 =    GetCovariance(4,1);
    const Double_t& mC32 =    GetCovariance(4,3);
    const Double_t& mC33 =    GetCovariance(4,4);
    const Double_t& mC44 = pV.GetCovariance(0,0);
    const Double_t& mC54 = pV.GetCovariance(1,0);
    const Double_t& mC55 = pV.GetCovariance(1,1);

    *timeErr2 =
      f5*mC55*f5 +
      f5*mC54*f4 +
      f4*mC44*f4 +
      f3*mC33*f3 +
      f3*mC32*f2 +
      f3*mC31*f1 +
      f3*mC30*f0 +
      f2*mC22*f2 +
      f2*mC21*f1 +
      f2*mC20*f0 +
      f1*mC11*f1 +
      f1*mC10*f0 +
      f0*mC00*f0;
  }
  return ( dx*Px() + dy*Py() )*mipt2;
}
 AliKFParticle.cxx:1
 AliKFParticle.cxx:2
 AliKFParticle.cxx:3
 AliKFParticle.cxx:4
 AliKFParticle.cxx:5
 AliKFParticle.cxx:6
 AliKFParticle.cxx:7
 AliKFParticle.cxx:8
 AliKFParticle.cxx:9
 AliKFParticle.cxx:10
 AliKFParticle.cxx:11
 AliKFParticle.cxx:12
 AliKFParticle.cxx:13
 AliKFParticle.cxx:14
 AliKFParticle.cxx:15
 AliKFParticle.cxx:16
 AliKFParticle.cxx:17
 AliKFParticle.cxx:18
 AliKFParticle.cxx:19
 AliKFParticle.cxx:20
 AliKFParticle.cxx:21
 AliKFParticle.cxx:22
 AliKFParticle.cxx:23
 AliKFParticle.cxx:24
 AliKFParticle.cxx:25
 AliKFParticle.cxx:26
 AliKFParticle.cxx:27
 AliKFParticle.cxx:28
 AliKFParticle.cxx:29
 AliKFParticle.cxx:30
 AliKFParticle.cxx:31
 AliKFParticle.cxx:32
 AliKFParticle.cxx:33
 AliKFParticle.cxx:34
 AliKFParticle.cxx:35
 AliKFParticle.cxx:36
 AliKFParticle.cxx:37
 AliKFParticle.cxx:38
 AliKFParticle.cxx:39
 AliKFParticle.cxx:40
 AliKFParticle.cxx:41
 AliKFParticle.cxx:42
 AliKFParticle.cxx:43
 AliKFParticle.cxx:44
 AliKFParticle.cxx:45
 AliKFParticle.cxx:46
 AliKFParticle.cxx:47
 AliKFParticle.cxx:48
 AliKFParticle.cxx:49
 AliKFParticle.cxx:50
 AliKFParticle.cxx:51
 AliKFParticle.cxx:52
 AliKFParticle.cxx:53
 AliKFParticle.cxx:54
 AliKFParticle.cxx:55
 AliKFParticle.cxx:56
 AliKFParticle.cxx:57
 AliKFParticle.cxx:58
 AliKFParticle.cxx:59
 AliKFParticle.cxx:60
 AliKFParticle.cxx:61
 AliKFParticle.cxx:62
 AliKFParticle.cxx:63
 AliKFParticle.cxx:64
 AliKFParticle.cxx:65
 AliKFParticle.cxx:66
 AliKFParticle.cxx:67
 AliKFParticle.cxx:68
 AliKFParticle.cxx:69
 AliKFParticle.cxx:70
 AliKFParticle.cxx:71
 AliKFParticle.cxx:72
 AliKFParticle.cxx:73
 AliKFParticle.cxx:74
 AliKFParticle.cxx:75
 AliKFParticle.cxx:76
 AliKFParticle.cxx:77
 AliKFParticle.cxx:78
 AliKFParticle.cxx:79
 AliKFParticle.cxx:80
 AliKFParticle.cxx:81
 AliKFParticle.cxx:82
 AliKFParticle.cxx:83
 AliKFParticle.cxx:84
 AliKFParticle.cxx:85
 AliKFParticle.cxx:86
 AliKFParticle.cxx:87
 AliKFParticle.cxx:88
 AliKFParticle.cxx:89
 AliKFParticle.cxx:90
 AliKFParticle.cxx:91
 AliKFParticle.cxx:92
 AliKFParticle.cxx:93
 AliKFParticle.cxx:94
 AliKFParticle.cxx:95
 AliKFParticle.cxx:96
 AliKFParticle.cxx:97
 AliKFParticle.cxx:98
 AliKFParticle.cxx:99
 AliKFParticle.cxx:100
 AliKFParticle.cxx:101
 AliKFParticle.cxx:102
 AliKFParticle.cxx:103
 AliKFParticle.cxx:104
 AliKFParticle.cxx:105
 AliKFParticle.cxx:106
 AliKFParticle.cxx:107
 AliKFParticle.cxx:108
 AliKFParticle.cxx:109
 AliKFParticle.cxx:110
 AliKFParticle.cxx:111
 AliKFParticle.cxx:112
 AliKFParticle.cxx:113
 AliKFParticle.cxx:114
 AliKFParticle.cxx:115
 AliKFParticle.cxx:116
 AliKFParticle.cxx:117
 AliKFParticle.cxx:118
 AliKFParticle.cxx:119
 AliKFParticle.cxx:120
 AliKFParticle.cxx:121
 AliKFParticle.cxx:122
 AliKFParticle.cxx:123
 AliKFParticle.cxx:124
 AliKFParticle.cxx:125
 AliKFParticle.cxx:126
 AliKFParticle.cxx:127
 AliKFParticle.cxx:128
 AliKFParticle.cxx:129
 AliKFParticle.cxx:130
 AliKFParticle.cxx:131
 AliKFParticle.cxx:132
 AliKFParticle.cxx:133
 AliKFParticle.cxx:134
 AliKFParticle.cxx:135
 AliKFParticle.cxx:136
 AliKFParticle.cxx:137
 AliKFParticle.cxx:138
 AliKFParticle.cxx:139
 AliKFParticle.cxx:140
 AliKFParticle.cxx:141
 AliKFParticle.cxx:142
 AliKFParticle.cxx:143
 AliKFParticle.cxx:144
 AliKFParticle.cxx:145
 AliKFParticle.cxx:146
 AliKFParticle.cxx:147
 AliKFParticle.cxx:148
 AliKFParticle.cxx:149
 AliKFParticle.cxx:150
 AliKFParticle.cxx:151
 AliKFParticle.cxx:152
 AliKFParticle.cxx:153
 AliKFParticle.cxx:154
 AliKFParticle.cxx:155
 AliKFParticle.cxx:156
 AliKFParticle.cxx:157
 AliKFParticle.cxx:158
 AliKFParticle.cxx:159
 AliKFParticle.cxx:160
 AliKFParticle.cxx:161
 AliKFParticle.cxx:162
 AliKFParticle.cxx:163
 AliKFParticle.cxx:164
 AliKFParticle.cxx:165
 AliKFParticle.cxx:166
 AliKFParticle.cxx:167
 AliKFParticle.cxx:168
 AliKFParticle.cxx:169
 AliKFParticle.cxx:170
 AliKFParticle.cxx:171
 AliKFParticle.cxx:172
 AliKFParticle.cxx:173
 AliKFParticle.cxx:174
 AliKFParticle.cxx:175
 AliKFParticle.cxx:176
 AliKFParticle.cxx:177
 AliKFParticle.cxx:178
 AliKFParticle.cxx:179
 AliKFParticle.cxx:180
 AliKFParticle.cxx:181
 AliKFParticle.cxx:182
 AliKFParticle.cxx:183
 AliKFParticle.cxx:184
 AliKFParticle.cxx:185
 AliKFParticle.cxx:186
 AliKFParticle.cxx:187
 AliKFParticle.cxx:188
 AliKFParticle.cxx:189
 AliKFParticle.cxx:190
 AliKFParticle.cxx:191
 AliKFParticle.cxx:192
 AliKFParticle.cxx:193
 AliKFParticle.cxx:194
 AliKFParticle.cxx:195
 AliKFParticle.cxx:196
 AliKFParticle.cxx:197
 AliKFParticle.cxx:198
 AliKFParticle.cxx:199
 AliKFParticle.cxx:200
 AliKFParticle.cxx:201
 AliKFParticle.cxx:202
 AliKFParticle.cxx:203
 AliKFParticle.cxx:204
 AliKFParticle.cxx:205
 AliKFParticle.cxx:206
 AliKFParticle.cxx:207
 AliKFParticle.cxx:208
 AliKFParticle.cxx:209
 AliKFParticle.cxx:210
 AliKFParticle.cxx:211
 AliKFParticle.cxx:212
 AliKFParticle.cxx:213
 AliKFParticle.cxx:214
 AliKFParticle.cxx:215
 AliKFParticle.cxx:216
 AliKFParticle.cxx:217
 AliKFParticle.cxx:218
 AliKFParticle.cxx:219
 AliKFParticle.cxx:220
 AliKFParticle.cxx:221
 AliKFParticle.cxx:222
 AliKFParticle.cxx:223
 AliKFParticle.cxx:224
 AliKFParticle.cxx:225
 AliKFParticle.cxx:226
 AliKFParticle.cxx:227
 AliKFParticle.cxx:228
 AliKFParticle.cxx:229
 AliKFParticle.cxx:230
 AliKFParticle.cxx:231
 AliKFParticle.cxx:232
 AliKFParticle.cxx:233
 AliKFParticle.cxx:234
 AliKFParticle.cxx:235
 AliKFParticle.cxx:236
 AliKFParticle.cxx:237
 AliKFParticle.cxx:238
 AliKFParticle.cxx:239
 AliKFParticle.cxx:240
 AliKFParticle.cxx:241
 AliKFParticle.cxx:242
 AliKFParticle.cxx:243
 AliKFParticle.cxx:244
 AliKFParticle.cxx:245
 AliKFParticle.cxx:246
 AliKFParticle.cxx:247
 AliKFParticle.cxx:248
 AliKFParticle.cxx:249
 AliKFParticle.cxx:250
 AliKFParticle.cxx:251
 AliKFParticle.cxx:252
 AliKFParticle.cxx:253
 AliKFParticle.cxx:254
 AliKFParticle.cxx:255
 AliKFParticle.cxx:256
 AliKFParticle.cxx:257
 AliKFParticle.cxx:258
 AliKFParticle.cxx:259
 AliKFParticle.cxx:260
 AliKFParticle.cxx:261
 AliKFParticle.cxx:262
 AliKFParticle.cxx:263
 AliKFParticle.cxx:264
 AliKFParticle.cxx:265
 AliKFParticle.cxx:266
 AliKFParticle.cxx:267
 AliKFParticle.cxx:268
 AliKFParticle.cxx:269
 AliKFParticle.cxx:270
 AliKFParticle.cxx:271
 AliKFParticle.cxx:272
 AliKFParticle.cxx:273
 AliKFParticle.cxx:274
 AliKFParticle.cxx:275
 AliKFParticle.cxx:276
 AliKFParticle.cxx:277
 AliKFParticle.cxx:278
 AliKFParticle.cxx:279
 AliKFParticle.cxx:280
 AliKFParticle.cxx:281
 AliKFParticle.cxx:282
 AliKFParticle.cxx:283
 AliKFParticle.cxx:284
 AliKFParticle.cxx:285
 AliKFParticle.cxx:286
 AliKFParticle.cxx:287
 AliKFParticle.cxx:288
 AliKFParticle.cxx:289
 AliKFParticle.cxx:290
 AliKFParticle.cxx:291
 AliKFParticle.cxx:292
 AliKFParticle.cxx:293
 AliKFParticle.cxx:294
 AliKFParticle.cxx:295
 AliKFParticle.cxx:296
 AliKFParticle.cxx:297
 AliKFParticle.cxx:298
 AliKFParticle.cxx:299
 AliKFParticle.cxx:300
 AliKFParticle.cxx:301
 AliKFParticle.cxx:302
 AliKFParticle.cxx:303
 AliKFParticle.cxx:304
 AliKFParticle.cxx:305
 AliKFParticle.cxx:306
 AliKFParticle.cxx:307
 AliKFParticle.cxx:308
 AliKFParticle.cxx:309
 AliKFParticle.cxx:310
 AliKFParticle.cxx:311
 AliKFParticle.cxx:312
 AliKFParticle.cxx:313
 AliKFParticle.cxx:314
 AliKFParticle.cxx:315
 AliKFParticle.cxx:316
 AliKFParticle.cxx:317
 AliKFParticle.cxx:318
 AliKFParticle.cxx:319
 AliKFParticle.cxx:320
 AliKFParticle.cxx:321
 AliKFParticle.cxx:322
 AliKFParticle.cxx:323
 AliKFParticle.cxx:324
 AliKFParticle.cxx:325
 AliKFParticle.cxx:326
 AliKFParticle.cxx:327
 AliKFParticle.cxx:328
 AliKFParticle.cxx:329
 AliKFParticle.cxx:330
 AliKFParticle.cxx:331
 AliKFParticle.cxx:332
 AliKFParticle.cxx:333
 AliKFParticle.cxx:334
 AliKFParticle.cxx:335
 AliKFParticle.cxx:336
 AliKFParticle.cxx:337
 AliKFParticle.cxx:338
 AliKFParticle.cxx:339
 AliKFParticle.cxx:340
 AliKFParticle.cxx:341
 AliKFParticle.cxx:342
 AliKFParticle.cxx:343
 AliKFParticle.cxx:344
 AliKFParticle.cxx:345
 AliKFParticle.cxx:346
 AliKFParticle.cxx:347
 AliKFParticle.cxx:348
 AliKFParticle.cxx:349
 AliKFParticle.cxx:350
 AliKFParticle.cxx:351
 AliKFParticle.cxx:352
 AliKFParticle.cxx:353
 AliKFParticle.cxx:354
 AliKFParticle.cxx:355
 AliKFParticle.cxx:356
 AliKFParticle.cxx:357
 AliKFParticle.cxx:358
 AliKFParticle.cxx:359
 AliKFParticle.cxx:360
 AliKFParticle.cxx:361
 AliKFParticle.cxx:362
 AliKFParticle.cxx:363
 AliKFParticle.cxx:364
 AliKFParticle.cxx:365
 AliKFParticle.cxx:366
 AliKFParticle.cxx:367
 AliKFParticle.cxx:368
 AliKFParticle.cxx:369
 AliKFParticle.cxx:370
 AliKFParticle.cxx:371
 AliKFParticle.cxx:372
 AliKFParticle.cxx:373
 AliKFParticle.cxx:374
 AliKFParticle.cxx:375
 AliKFParticle.cxx:376
 AliKFParticle.cxx:377
 AliKFParticle.cxx:378
 AliKFParticle.cxx:379
 AliKFParticle.cxx:380
 AliKFParticle.cxx:381
 AliKFParticle.cxx:382
 AliKFParticle.cxx:383
 AliKFParticle.cxx:384
 AliKFParticle.cxx:385
 AliKFParticle.cxx:386
 AliKFParticle.cxx:387
 AliKFParticle.cxx:388
 AliKFParticle.cxx:389
 AliKFParticle.cxx:390
 AliKFParticle.cxx:391
 AliKFParticle.cxx:392
 AliKFParticle.cxx:393
 AliKFParticle.cxx:394
 AliKFParticle.cxx:395
 AliKFParticle.cxx:396
 AliKFParticle.cxx:397
 AliKFParticle.cxx:398
 AliKFParticle.cxx:399
 AliKFParticle.cxx:400
 AliKFParticle.cxx:401
 AliKFParticle.cxx:402
 AliKFParticle.cxx:403
 AliKFParticle.cxx:404
 AliKFParticle.cxx:405
 AliKFParticle.cxx:406
 AliKFParticle.cxx:407
 AliKFParticle.cxx:408
 AliKFParticle.cxx:409
 AliKFParticle.cxx:410
 AliKFParticle.cxx:411
 AliKFParticle.cxx:412
 AliKFParticle.cxx:413
 AliKFParticle.cxx:414
 AliKFParticle.cxx:415
 AliKFParticle.cxx:416
 AliKFParticle.cxx:417
 AliKFParticle.cxx:418
 AliKFParticle.cxx:419
 AliKFParticle.cxx:420
 AliKFParticle.cxx:421
 AliKFParticle.cxx:422
 AliKFParticle.cxx:423
 AliKFParticle.cxx:424
 AliKFParticle.cxx:425
 AliKFParticle.cxx:426
 AliKFParticle.cxx:427
 AliKFParticle.cxx:428
 AliKFParticle.cxx:429
 AliKFParticle.cxx:430
 AliKFParticle.cxx:431
 AliKFParticle.cxx:432
 AliKFParticle.cxx:433
 AliKFParticle.cxx:434
 AliKFParticle.cxx:435
 AliKFParticle.cxx:436
 AliKFParticle.cxx:437
 AliKFParticle.cxx:438
 AliKFParticle.cxx:439
 AliKFParticle.cxx:440
 AliKFParticle.cxx:441
 AliKFParticle.cxx:442
 AliKFParticle.cxx:443
 AliKFParticle.cxx:444
 AliKFParticle.cxx:445
 AliKFParticle.cxx:446
 AliKFParticle.cxx:447
 AliKFParticle.cxx:448
 AliKFParticle.cxx:449
 AliKFParticle.cxx:450
 AliKFParticle.cxx:451
 AliKFParticle.cxx:452
 AliKFParticle.cxx:453
 AliKFParticle.cxx:454
 AliKFParticle.cxx:455
 AliKFParticle.cxx:456
 AliKFParticle.cxx:457
 AliKFParticle.cxx:458
 AliKFParticle.cxx:459
 AliKFParticle.cxx:460
 AliKFParticle.cxx:461
 AliKFParticle.cxx:462
 AliKFParticle.cxx:463
 AliKFParticle.cxx:464
 AliKFParticle.cxx:465
 AliKFParticle.cxx:466
 AliKFParticle.cxx:467
 AliKFParticle.cxx:468
 AliKFParticle.cxx:469
 AliKFParticle.cxx:470
 AliKFParticle.cxx:471
 AliKFParticle.cxx:472
 AliKFParticle.cxx:473
 AliKFParticle.cxx:474
 AliKFParticle.cxx:475
 AliKFParticle.cxx:476
 AliKFParticle.cxx:477
 AliKFParticle.cxx:478
 AliKFParticle.cxx:479
 AliKFParticle.cxx:480
 AliKFParticle.cxx:481
 AliKFParticle.cxx:482
 AliKFParticle.cxx:483
 AliKFParticle.cxx:484
 AliKFParticle.cxx:485
 AliKFParticle.cxx:486
 AliKFParticle.cxx:487
 AliKFParticle.cxx:488
 AliKFParticle.cxx:489
 AliKFParticle.cxx:490
 AliKFParticle.cxx:491