]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliKFParticleBase.cxx
Compilation on Windoiws/Cygwin
[u/mrichter/AliRoot.git] / STEER / AliKFParticleBase.cxx
1 //---------------------------------------------------------------------------------
2 // Implementation of the AliKFParticleBase class
3 // .
4 // @author  S.Gorbunov, I.Kisel
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 describes general mathematics which is used by AliKFParticle class
14 // 
15 //  -= Copyright &copy ALICE HLT Group =-
16 //_________________________________________________________________________________
17
18
19 #include "AliKFParticleBase.h"
20 #include "TMath.h"
21
22 ClassImp(AliKFParticleBase)
23
24
25 AliKFParticleBase::AliKFParticleBase() :fQ(0), fNDF(-3), fChi2(0), fSFromDecay(0), fAtProductionVertex(0), fIsLinearized(0)
26
27   //* Constructor 
28
29   Initialize();
30 }
31
32 void AliKFParticleBase::Initialize( const Double_t Param[], const Double_t Cov[], Int_t Charge, Double_t Mass )
33 {
34   // Constructor from "cartesian" track, particle mass hypothesis should be provided
35   //
36   // Param[6] = { X, Y, Z, Px, Py, Pz } - position and momentum
37   // Cov [21] = lower-triangular part of the covariance matrix:
38   //
39   //                (  0  .  .  .  .  . )
40   //                (  1  2  .  .  .  . )
41   //  Cov. matrix = (  3  4  5  .  .  . ) - numbering of covariance elements in Cov[]
42   //                (  6  7  8  9  .  . )
43   //                ( 10 11 12 13 14  . )
44   //                ( 15 16 17 18 19 20 )
45
46
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];
49
50   Double_t energy = TMath::Sqrt( Mass*Mass + fP[3]*fP[3] + fP[4]*fP[4] + fP[5]*fP[5]);
51   fP[6] = energy;
52   fP[7] = 0;
53   fQ = Charge;
54   fNDF = 0;
55   fChi2 = 0;
56   fAtProductionVertex = 0;
57   fIsLinearized = 0;
58   fSFromDecay = 0;
59
60   Double_t energyInv = 1./energy;
61   Double_t 
62     h0 = fP[3]*energyInv,
63     h1 = fP[4]*energyInv,
64     h2 = fP[5]*energyInv;
65
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;
75   fC[35] = 1.;
76 }
77
78 void AliKFParticleBase::Initialize()
79 {
80   //* Initialise covariance matrix and set current parameters to 0.0 
81
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.;
85   fC[35] = 1.;
86   fNDF  = -3;
87   fChi2 =  0.;
88   fQ = 0;
89   fSFromDecay = 0;
90   fAtProductionVertex = 0;
91   fVtxGuess[0]=fVtxGuess[1]=fVtxGuess[2]=0.;
92   fIsLinearized = 0;
93 }
94
95 void AliKFParticleBase::SetVtxGuess( Double_t x, Double_t y, Double_t z )
96 {
97   //* Set decay vertex parameters for linearisation 
98
99   fVtxGuess[0] = x;
100   fVtxGuess[1] = y;
101   fVtxGuess[2] = z;
102   fIsLinearized = 1;
103 }
104
105
106 Int_t AliKFParticleBase::GetMomentum( Double_t &pmom, Double_t &sigmap )  const 
107 {
108   //* Calculate particle momentum
109
110   Double_t x = fP[3];
111   Double_t y = fP[4];
112   Double_t z = fP[5];
113   Double_t x2 = x*x;
114   Double_t y2 = y*y;
115   Double_t z2 = z*z;
116   Double_t p2 = x2+y2+z2;
117   pmom = TMath::Sqrt(p2);
118   sigmap = (x2*fC[9]+y2*fC[14]+z2*fC[20] + 2*(x*y*fC[13]+x*z*fC[18]+y*z*fC[19]) );
119   if( sigmap>0 && pmom>1.e-4 ){
120     sigmap = TMath::Sqrt(sigmap)/pmom;
121     return 0;
122   }
123   return 1;
124 }
125
126 Int_t AliKFParticleBase::GetMass( Double_t &mass, Double_t &sigmam ) const 
127 {
128   //* Calculate particle mass
129   
130   // s = sigma^2 of m2/2
131
132   Double_t s = (  fP[3]*fP[3]*fC[9] + fP[4]*fP[4]*fC[14] + fP[5]*fP[5]*fC[20] 
133                   + fP[6]*fP[6]*fC[27] 
134                 +2*( + fP[3]*fP[4]*fC[13] + fP[5]*(fP[3]*fC[18] + fP[4]*fC[19]) 
135                      - fP[6]*( fP[3]*fC[24] + fP[4]*fC[25] + fP[5]*fC[26] )   )
136                  ); 
137   Double_t m2 = fP[6]*fP[6] - fP[3]*fP[3] - fP[4]*fP[4] - fP[5]*fP[5];
138   mass     = 0;
139   if( m2>1.e-20 ){
140     mass = TMath::Sqrt(m2);
141     if( s>=0 ){
142       sigmam = TMath::Sqrt(s/m2);
143       return 0;
144     }
145   }
146   sigmam = 1.e20;
147   return 1;
148 }
149
150
151 Int_t AliKFParticleBase::GetDecayLength( Double_t &dlen, Double_t &sigmal ) const 
152 {
153   //* Calculate particle decay length [cm]
154
155   Double_t x = fP[3];
156   Double_t y = fP[4];
157   Double_t z = fP[5];
158   Double_t t = fP[7];
159   Double_t x2 = x*x;
160   Double_t y2 = y*y;
161   Double_t z2 = z*z;
162   Double_t p2 = x2+y2+z2;
163   dlen = t*TMath::Sqrt(p2);
164   if( p2>1.e-4){
165     sigmal = p2*fC[35] + t*t/p2*(x2*fC[9]+y2*fC[14]+z2*fC[20]
166                                 + 2*(x*y*fC[13]+x*z*fC[18]+y*z*fC[19]) )
167       + 2*t*(x*fC[31]+y*fC[32]+z*fC[33]);
168     sigmal = TMath::Sqrt(TMath::Abs(sigmal));
169     return 0;
170   }
171   sigmal = 1.e20;
172   return 1;
173 }
174
175 Int_t AliKFParticleBase::GetLifeTime( Double_t &tauc, Double_t &sigmat ) const 
176 {
177   //* Calculate particle decay time [s]
178
179   Double_t m, dm;
180   GetMass( m, dm );
181   Double_t cTM = (-fP[3]*fC[31] - fP[4]*fC[32] - fP[5]*fC[33] + fP[6]*fC[34]);
182   tauc = fP[7]*m;
183   sigmat = m*m*fC[35] + 2*fP[7]*cTM + fP[7]*fP[7]*dm*dm;
184   if( sigmat > 0 ){
185     sigmat = TMath::Sqrt( sigmat );
186     return 0;
187   }
188   sigmat = 1.e20;
189   return 1;
190 }
191
192
193 void AliKFParticleBase::operator +=( const AliKFParticleBase &Daughter )
194 {
195   //* Add daughter via operator+=
196
197   AddDaughter( Daughter );
198 }
199   
200
201 void AliKFParticleBase::GetMeasurement( const Double_t XYZ[], Double_t m[], Double_t V[] ) const
202 {
203   //* Get additional covariances V used during measurement
204
205   Double_t b[3];
206   GetFieldValue( XYZ, b );
207   const Double_t kCLight =  0.000299792458;
208   b[0]*=kCLight; b[1]*=kCLight; b[2]*=kCLight;
209
210   Transport( GetDStoPoint(XYZ), m, V );
211
212   Double_t d[3] = { XYZ[0]-m[0], XYZ[1]-m[1], XYZ[2]-m[2] };
213   Double_t p2 = m[3]*m[3]+m[4]*m[4]+m[5]*m[5];
214   Double_t sigmaS = (p2>1.e-4) ? ( .1+TMath::Sqrt( 100*(d[0]*d[0]+d[1]*d[1]+d[2]*d[2])) )/TMath::Sqrt(p2) : 1.;
215
216   Double_t h[6];
217
218   h[0] = m[3]*sigmaS;
219   h[1] = m[4]*sigmaS;
220   h[2] = m[5]*sigmaS;
221   h[3] = ( h[1]*b[2]-h[2]*b[1] )*GetQ();
222   h[4] = ( h[2]*b[0]-h[0]*b[2] )*GetQ();
223   h[5] = ( h[0]*b[1]-h[1]*b[0] )*GetQ();
224     
225   //* Fit of momentum (Px,Py,Pz) to XYZ point
226   if(0){      
227     Double_t mVv[6] = 
228       { V[ 0] + h[0]*h[0], 
229         V[ 1] + h[0]*h[1], V[ 2] + h[1]*h[1],                        
230         V[ 3] + h[0]*h[2], V[ 4] + h[1]*h[2], V[ 5] + h[2]*h[2] };
231     
232     Double_t mVvp[9]=
233       { V[ 6] + h[0]*h[3], V[ 7] + h[1]*h[3], V[ 8] + h[2]*h[3],
234         V[10] + h[0]*h[4], V[11] + h[1]*h[4], V[12] + h[2]*h[4],
235         V[15] + h[0]*h[5], V[16] + h[1]*h[5], V[17] + h[2]*h[5] };
236     
237     Double_t mS[6] = 
238       { mVv[2]*mVv[5] - mVv[4]*mVv[4],
239         mVv[3]*mVv[4] - mVv[1]*mVv[5], mVv[0]*mVv[5] - mVv[3]*mVv[3],
240         mVv[1]*mVv[4] - mVv[2]*mVv[3], mVv[1]*mVv[3] - mVv[0]*mVv[4], 
241         mVv[0]*mVv[2] - mVv[1]*mVv[1] };                 
242     
243     Double_t s = ( mVv[0]*mS[0] + mVv[1]*mS[1] + mVv[3]*mS[3] );    
244     s = ( s > 1.E-20 )  ?1./s :0;
245         
246     Double_t mSz[3] = { mS[0]*d[0]+mS[1]*d[1]+mS[3]*d[2],
247                         mS[1]*d[0]+mS[2]*d[1]+mS[4]*d[2],
248                         mS[3]*d[0]+mS[4]*d[1]+mS[5]*d[2] };
249     
250     Double_t px = m[3] + s*( mVvp[0]*mSz[0] + mVvp[1]*mSz[1] + mVvp[2]*mSz[2] );
251     Double_t py = m[4] + s*( mVvp[3]*mSz[0] + mVvp[4]*mSz[1] + mVvp[5]*mSz[2] );
252     Double_t pz = m[5] + s*( mVvp[6]*mSz[0] + mVvp[7]*mSz[1] + mVvp[8]*mSz[2] );
253
254     h[0] = px*sigmaS;
255     h[1] = py*sigmaS;
256     h[2] = pz*sigmaS;
257     h[3] = ( h[1]*b[2]-h[2]*b[1] )*GetQ();
258     h[4] = ( h[2]*b[0]-h[0]*b[2] )*GetQ();
259     h[5] = ( h[0]*b[1]-h[1]*b[0] )*GetQ();
260   }
261    
262   V[ 0]+= h[0]*h[0];
263   V[ 1]+= h[1]*h[0];
264   V[ 2]+= h[1]*h[1];
265   V[ 3]+= h[2]*h[0];
266   V[ 4]+= h[2]*h[1];
267   V[ 5]+= h[2]*h[2];
268
269   V[ 6]+= h[3]*h[0];
270   V[ 7]+= h[3]*h[1];
271   V[ 8]+= h[3]*h[2];
272   V[ 9]+= h[3]*h[3];
273
274   V[10]+= h[4]*h[0];
275   V[11]+= h[4]*h[1];
276   V[12]+= h[4]*h[2];
277   V[13]+= h[4]*h[3];
278   V[14]+= h[4]*h[4];
279
280   V[15]+= h[5]*h[0];
281   V[16]+= h[5]*h[1];
282   V[17]+= h[5]*h[2];
283   V[18]+= h[5]*h[3];
284   V[19]+= h[5]*h[4];
285   V[20]+= h[5]*h[5];
286 }
287
288   
289 void AliKFParticleBase::AddDaughter( const AliKFParticleBase &Daughter )
290 {
291   //* Add daughter 
292
293   if( fNDF<-1 ){ // first daughter -> just copy
294     fNDF   = -1;
295     fQ     =  Daughter.GetQ();
296     for( Int_t i=0; i<7; i++) fP[i] = Daughter.fP[i];
297     for( Int_t i=0; i<28; i++) fC[i] = Daughter.fC[i];
298     fSFromDecay = 0;
299     return;
300   }
301
302   TransportToDecayVertex();
303
304   Double_t b[3]; 
305   Int_t maxIter = 1;
306
307   if( !fIsLinearized ){
308     if( fNDF==-1 ){
309       Double_t ds, ds1;
310       GetDStoParticle(Daughter, ds, ds1);      
311       TransportToDS( ds );
312       Double_t m[8];
313       Double_t mCd[36];       
314       Daughter.Transport( ds1, m, mCd );    
315       fVtxGuess[0] = .5*( fP[0] + m[0] );
316       fVtxGuess[1] = .5*( fP[1] + m[1] );
317       fVtxGuess[2] = .5*( fP[2] + m[2] );
318     } else {
319       fVtxGuess[0] = fP[0];
320       fVtxGuess[1] = fP[1];
321       fVtxGuess[2] = fP[2]; 
322     }
323     maxIter = 3;
324   }
325
326   for( Int_t iter=0; iter<maxIter; iter++ ){
327
328     {
329       GetFieldValue( fVtxGuess, b );
330       const Double_t kCLight =  0.000299792458;
331       b[0]*=kCLight; b[1]*=kCLight; b[2]*=kCLight;
332     }
333
334     Double_t *ffP = fP, *ffC = fC, tmpP[8], tmpC[36];
335     if( fNDF==-1 ){      
336       GetMeasurement( fVtxGuess, tmpP, tmpC );
337       ffP = tmpP;
338       ffC = tmpC;
339     }
340
341     Double_t m[8], mV[36];
342     if( Daughter.fC[35]>0 ){
343       Daughter.GetMeasurement( fVtxGuess, m, mV );
344     } else {
345       for( Int_t i=0; i<8; i++ ) m[i] = Daughter.fP[i];
346       for( Int_t i=0; i<36; i++ ) mV[i] = Daughter.fC[i];
347     }
348
349     //*
350
351     Double_t mS[6];
352     {
353       Double_t mSi[6] = { ffC[0]+mV[0], 
354                           ffC[1]+mV[1], ffC[2]+mV[2], 
355                           ffC[3]+mV[3], ffC[4]+mV[4], ffC[5]+mV[5] };
356       
357       mS[0] = mSi[2]*mSi[5] - mSi[4]*mSi[4];
358       mS[1] = mSi[3]*mSi[4] - mSi[1]*mSi[5];
359       mS[2] = mSi[0]*mSi[5] - mSi[3]*mSi[3];
360       mS[3] = mSi[1]*mSi[4] - mSi[2]*mSi[3];
361       mS[4] = mSi[1]*mSi[3] - mSi[0]*mSi[4];
362       mS[5] = mSi[0]*mSi[2] - mSi[1]*mSi[1];     
363       
364       Double_t s = ( mSi[0]*mS[0] + mSi[1]*mS[1] + mSi[3]*mS[3] );
365       s = ( s > 1.E-20 )  ?1./s :0;       
366       mS[0]*=s;
367       mS[1]*=s;
368       mS[2]*=s;
369       mS[3]*=s;
370       mS[4]*=s;
371       mS[5]*=s;
372     }
373     
374     //* Residual (measured - estimated)
375     
376     Double_t zeta[3] = { m[0]-ffP[0], m[1]-ffP[1], m[2]-ffP[2] };    
377     
378     //* CHt = CH' - D'
379     
380     Double_t mCHt0[7], mCHt1[7], mCHt2[7];
381     
382     mCHt0[0]=ffC[ 0] ;       mCHt1[0]=ffC[ 1] ;       mCHt2[0]=ffC[ 3] ;
383     mCHt0[1]=ffC[ 1] ;       mCHt1[1]=ffC[ 2] ;       mCHt2[1]=ffC[ 4] ;
384     mCHt0[2]=ffC[ 3] ;       mCHt1[2]=ffC[ 4] ;       mCHt2[2]=ffC[ 5] ;
385     mCHt0[3]=ffC[ 6]-mV[ 6]; mCHt1[3]=ffC[ 7]-mV[ 7]; mCHt2[3]=ffC[ 8]-mV[ 8];
386     mCHt0[4]=ffC[10]-mV[10]; mCHt1[4]=ffC[11]-mV[11]; mCHt2[4]=ffC[12]-mV[12];
387     mCHt0[5]=ffC[15]-mV[15]; mCHt1[5]=ffC[16]-mV[16]; mCHt2[5]=ffC[17]-mV[17];
388     mCHt0[6]=ffC[21]-mV[21]; mCHt1[6]=ffC[22]-mV[22]; mCHt2[6]=ffC[23]-mV[23];
389   
390     //* Kalman gain K = mCH'*S
391     
392     Double_t k0[7], k1[7], k2[7];
393     
394     for(Int_t i=0;i<7;++i){
395       k0[i] = mCHt0[i]*mS[0] + mCHt1[i]*mS[1] + mCHt2[i]*mS[3];
396       k1[i] = mCHt0[i]*mS[1] + mCHt1[i]*mS[2] + mCHt2[i]*mS[4];
397       k2[i] = mCHt0[i]*mS[3] + mCHt1[i]*mS[4] + mCHt2[i]*mS[5];
398     }
399
400    //* New estimation of the vertex position 
401
402     if( iter<maxIter-1 ){
403       for(Int_t i=0; i<3; ++i) 
404         fVtxGuess[i]= ffP[i] + k0[i]*zeta[0]+k1[i]*zeta[1]+k2[i]*zeta[2];
405       continue;
406     }
407
408     // last itearation -> update the particle
409
410     //* Add the daughter momentum to the particle momentum
411     
412     ffP[ 3] += m[ 3];
413     ffP[ 4] += m[ 4];
414     ffP[ 5] += m[ 5];
415     ffP[ 6] += m[ 6];
416   
417     ffC[ 9] += mV[ 9];
418     ffC[13] += mV[13];
419     ffC[14] += mV[14];
420     ffC[18] += mV[18];
421     ffC[19] += mV[19];
422     ffC[20] += mV[20];
423     ffC[24] += mV[24];
424     ffC[25] += mV[25];
425     ffC[26] += mV[26];
426     ffC[27] += mV[27];
427     
428     //* New estimation of the vertex position r += K*zeta
429     
430     for(Int_t i=0;i<7;++i) 
431       fP[i] = ffP[i] + k0[i]*zeta[0] + k1[i]*zeta[1] + k2[i]*zeta[2];
432     
433     //* New covariance matrix C -= K*(mCH')'
434    
435     for(Int_t i=0, k=0;i<7;++i){
436       for(Int_t j=0;j<=i;++j,++k) 
437         fC[k] = ffC[k] - (k0[i]*mCHt0[j] + k1[i]*mCHt1[j] + k2[i]*mCHt2[j] );
438     }
439     
440     //* Calculate Chi^2 
441
442     fNDF  += 2;
443     fQ    +=  Daughter.GetQ();
444     fSFromDecay = 0;    
445     fChi2 += (mS[0]*zeta[0] + mS[1]*zeta[1] + mS[3]*zeta[2])*zeta[0]
446       +      (mS[1]*zeta[0] + mS[2]*zeta[1] + mS[4]*zeta[2])*zeta[1]
447       +      (mS[3]*zeta[0] + mS[4]*zeta[1] + mS[5]*zeta[2])*zeta[2];  
448   }
449 }
450
451 void AliKFParticleBase::SetProductionVertex( const AliKFParticleBase &Vtx )
452 {
453   //* Set production vertex for the particle
454
455   const Double_t *m = Vtx.fP, *mV = Vtx.fC;
456
457   Bool_t noS = ( fC[35]<=0 ); // no decay length allowed
458
459   if( noS ){ 
460     TransportToDecayVertex();
461     fP[7] = 0;
462     fC[28] = fC[29] = fC[30] = fC[31] = fC[32] = fC[33] = fC[35] = fC[35] = 0;
463   } else {
464     TransportToDS( GetDStoPoint( m ) );
465     fP[7] = -fSFromDecay;
466     Convert(1);
467   }
468
469  
470   Double_t mAi[6];
471   mAi[0] = fC[4]*fC[4] - fC[2]*fC[5];
472   mAi[1] = fC[1]*fC[5] - fC[3]*fC[4];
473   mAi[3] = fC[2]*fC[3] - fC[1]*fC[4];
474   Double_t det = 1./(fC[0]*mAi[0] + fC[1]*mAi[1] + fC[3]*mAi[3]);
475   mAi[0] *= det;
476   mAi[1] *= det;
477   mAi[3] *= det;
478   mAi[2] = ( fC[3]*fC[3] - fC[0]*fC[5] )*det;
479   mAi[4] = ( fC[0]*fC[4] - fC[1]*fC[3] )*det;
480   mAi[5] = ( fC[1]*fC[1] - fC[0]*fC[2] )*det;
481
482   Double_t mB[5][3];
483
484   mB[0][0] = fC[ 6]*mAi[0] + fC[ 7]*mAi[1] + fC[ 8]*mAi[3];
485   mB[0][1] = fC[ 6]*mAi[1] + fC[ 7]*mAi[2] + fC[ 8]*mAi[4];
486   mB[0][2] = fC[ 6]*mAi[3] + fC[ 7]*mAi[4] + fC[ 8]*mAi[5];
487
488   mB[1][0] = fC[10]*mAi[0] + fC[11]*mAi[1] + fC[12]*mAi[3];
489   mB[1][1] = fC[10]*mAi[1] + fC[11]*mAi[2] + fC[12]*mAi[4];
490   mB[1][2] = fC[10]*mAi[3] + fC[11]*mAi[4] + fC[12]*mAi[5];
491
492   mB[2][0] = fC[15]*mAi[0] + fC[16]*mAi[1] + fC[17]*mAi[3];
493   mB[2][1] = fC[15]*mAi[1] + fC[16]*mAi[2] + fC[17]*mAi[4];
494   mB[2][2] = fC[15]*mAi[3] + fC[16]*mAi[4] + fC[17]*mAi[5];
495
496   mB[3][0] = fC[21]*mAi[0] + fC[22]*mAi[1] + fC[23]*mAi[3];
497   mB[3][1] = fC[21]*mAi[1] + fC[22]*mAi[2] + fC[23]*mAi[4];
498   mB[3][2] = fC[21]*mAi[3] + fC[22]*mAi[4] + fC[23]*mAi[5];
499
500   mB[4][0] = fC[28]*mAi[0] + fC[29]*mAi[1] + fC[30]*mAi[3];
501   mB[4][1] = fC[28]*mAi[1] + fC[29]*mAi[2] + fC[30]*mAi[4];
502   mB[4][2] = fC[28]*mAi[3] + fC[29]*mAi[4] + fC[30]*mAi[5];
503
504   Double_t z[3] = { m[0]-fP[0], m[1]-fP[1], m[2]-fP[2] };
505
506   {
507     Double_t mAV[6] = { fC[0]-mV[0], fC[1]-mV[1], fC[2]-mV[2], 
508                         fC[3]-mV[3], fC[4]-mV[4], fC[5]-mV[5] };
509     Double_t mAVi[6];
510     mAVi[0] = mAV[4]*mAV[4] - mAV[2]*mAV[5];
511     mAVi[1] = mAV[1]*mAV[5] - mAV[3]*mAV[4];
512     mAVi[2] = mAV[3]*mAV[3] - mAV[0]*mAV[5] ;
513     mAVi[3] = mAV[2]*mAV[3] - mAV[1]*mAV[4];
514     mAVi[4] = mAV[0]*mAV[4] - mAV[1]*mAV[3] ;
515     mAVi[5] = mAV[1]*mAV[1] - mAV[0]*mAV[2] ;
516
517     det = ( mAV[0]*mAVi[0] + mAV[1]*mAVi[1] + mAV[3]*mAVi[3] );
518     if( det>1.e-8 ) det = 1./det;
519
520     fNDF  += 2;
521     fChi2 += 
522       ( +(mAVi[0]*z[0] + mAVi[1]*z[1] + mAVi[3]*z[2])*z[0]
523         +(mAVi[1]*z[0] + mAVi[2]*z[1] + mAVi[4]*z[2])*z[1]
524         +(mAVi[3]*z[0] + mAVi[4]*z[1] + mAVi[5]*z[2])*z[2] )*det;
525   }
526   
527   fP[0] = m[0];
528   fP[1] = m[1];
529   fP[2] = m[2];
530   fP[3]+= mB[0][0]*z[0] + mB[0][1]*z[1] + mB[0][2]*z[2];
531   fP[4]+= mB[1][0]*z[0] + mB[1][1]*z[1] + mB[1][2]*z[2];
532   fP[5]+= mB[2][0]*z[0] + mB[2][1]*z[1] + mB[2][2]*z[2];
533   fP[6]+= mB[3][0]*z[0] + mB[3][1]*z[1] + mB[3][2]*z[2];
534   fP[7]+= mB[4][0]*z[0] + mB[4][1]*z[1] + mB[4][2]*z[2];
535   
536   Double_t d0, d1, d2;
537
538   fC[0] = mV[0];
539   fC[1] = mV[1];
540   fC[2] = mV[2];
541   fC[3] = mV[3];
542   fC[4] = mV[4];
543   fC[5] = mV[5];
544
545   d0= mB[0][0]*mV[0] + mB[0][1]*mV[1] + mB[0][2]*mV[3] - fC[ 6];
546   d1= mB[0][0]*mV[1] + mB[0][1]*mV[2] + mB[0][2]*mV[4] - fC[ 7];
547   d2= mB[0][0]*mV[3] + mB[0][1]*mV[4] + mB[0][2]*mV[5] - fC[ 8];
548
549   fC[ 6]+= d0;
550   fC[ 7]+= d1;
551   fC[ 8]+= d2;
552   fC[ 9]+= d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
553
554   d0= mB[1][0]*mV[0] + mB[1][1]*mV[1] + mB[1][2]*mV[3] - fC[10];
555   d1= mB[1][0]*mV[1] + mB[1][1]*mV[2] + mB[1][2]*mV[4] - fC[11];
556   d2= mB[1][0]*mV[3] + mB[1][1]*mV[4] + mB[1][2]*mV[5] - fC[12];
557
558   fC[10]+= d0;
559   fC[11]+= d1;
560   fC[12]+= d2;
561   fC[13]+= d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
562   fC[14]+= d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
563
564   d0= mB[2][0]*mV[0] + mB[2][1]*mV[1] + mB[2][2]*mV[3] - fC[15];
565   d1= mB[2][0]*mV[1] + mB[2][1]*mV[2] + mB[2][2]*mV[4] - fC[16];
566   d2= mB[2][0]*mV[3] + mB[2][1]*mV[4] + mB[2][2]*mV[5] - fC[17];
567
568   fC[15]+= d0;
569   fC[16]+= d1;
570   fC[17]+= d2;
571   fC[18]+= d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
572   fC[19]+= d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
573   fC[20]+= d0*mB[2][0] + d1*mB[2][1] + d2*mB[2][2];
574
575   d0= mB[3][0]*mV[0] + mB[3][1]*mV[1] + mB[3][2]*mV[3] - fC[21];
576   d1= mB[3][0]*mV[1] + mB[3][1]*mV[2] + mB[3][2]*mV[4] - fC[22];
577   d2= mB[3][0]*mV[3] + mB[3][1]*mV[4] + mB[3][2]*mV[5] - fC[23];
578
579   fC[21]+= d0;
580   fC[22]+= d1;
581   fC[23]+= d2;
582   fC[24]+= d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
583   fC[25]+= d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
584   fC[26]+= d0*mB[2][0] + d1*mB[2][1] + d2*mB[2][2];
585   fC[27]+= d0*mB[3][0] + d1*mB[3][1] + d2*mB[3][2];
586
587   d0= mB[4][0]*mV[0] + mB[4][1]*mV[1] + mB[4][2]*mV[3] - fC[28];
588   d1= mB[4][0]*mV[1] + mB[4][1]*mV[2] + mB[4][2]*mV[4] - fC[29];
589   d2= mB[4][0]*mV[3] + mB[4][1]*mV[4] + mB[4][2]*mV[5] - fC[30];
590
591   fC[28]+= d0;
592   fC[29]+= d1;
593   fC[30]+= d2;
594   fC[31]+= d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
595   fC[32]+= d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
596   fC[33]+= d0*mB[2][0] + d1*mB[2][1] + d2*mB[2][2];
597   fC[34]+= d0*mB[3][0] + d1*mB[3][1] + d2*mB[3][2];
598   fC[35]+= d0*mB[4][0] + d1*mB[4][1] + d2*mB[4][2];
599   
600   if( noS ){ 
601     fP[7] = 0;
602     fC[28] = fC[29] = fC[30] = fC[31] = fC[32] = fC[33] = fC[35] = fC[35] = 0;
603   } else {
604     TransportToDS( fP[7] );
605     Convert(0);
606   }
607
608   fSFromDecay = 0;
609 }
610
611
612
613
614 void AliKFParticleBase::SetMassConstraint( Double_t Mass, Double_t SigmaMass )
615 {  
616   //* Set hard mass constraint 
617
618   Double_t mH[8];
619   mH[0] = mH[1] = mH[2] = 0.;
620   mH[3] = -2*fP[3]; 
621   mH[4] = -2*fP[4]; 
622   mH[5] = -2*fP[5]; 
623   mH[6] =  2*fP[6];
624   mH[7] = 0; 
625   Double_t m2 = ( fP[6]*fP[6] 
626                 - fP[3]*fP[3] - fP[4]*fP[4] - fP[5]*fP[5] ); 
627
628   Double_t zeta = Mass*Mass - m2;
629   for(Int_t i=0;i<8;++i) zeta -= mH[i]*(fP[i]-fP[i]);
630   
631   Double_t s = 4*Mass*Mass*SigmaMass*SigmaMass;
632   Double_t mCHt[8];
633   for (Int_t i=0;i<8;++i ){
634     mCHt[i] = 0.0;
635     for (Int_t j=0;j<8;++j) mCHt[i] += Cij(i,j)*mH[j];
636     s += mH[i]*mCHt[i];
637   }
638   
639   if( s<1.e-20 ) return;
640   s = 1./s;
641   fChi2 += zeta*zeta*s;
642   fNDF  += 1;
643   for( Int_t i=0, ii=0; i<8; ++i ){
644     Double_t ki = mCHt[i]*s;
645     fP[i]+= ki*zeta;
646     for(Int_t j=0;j<=i;++j) fC[ii++] -= ki*mCHt[j];    
647   }
648 }
649
650
651 void AliKFParticleBase::SetNoDecayLength()
652 {  
653   //* Set no decay length for resonances
654
655   TransportToDecayVertex();
656
657   Double_t h[8];
658   h[0] = h[1] = h[2] = h[3] = h[4] = h[5] = h[6] = 0;
659   h[7] = 1; 
660
661   Double_t zeta = 0 - fP[7];
662   for(Int_t i=0;i<8;++i) zeta -= h[i]*(fP[i]-fP[i]);
663   
664   Double_t s = fC[35];   
665   if( s>1.e-20 ){
666     s = 1./s;
667     fChi2 += zeta*zeta*s;
668     fNDF  += 1;
669     for( Int_t i=0, ii=0; i<7; ++i ){
670       Double_t ki = fC[28+i]*s;
671       fP[i]+= ki*zeta;
672       for(Int_t j=0;j<=i;++j) fC[ii++] -= ki*fC[28+j];    
673     }
674   }
675   fP[7] = 0;
676   fC[28] = fC[29] = fC[30] = fC[31] = fC[32] = fC[33] = fC[35] = fC[35] = 0;
677 }
678
679
680 void AliKFParticleBase::Construct( const AliKFParticleBase* vDaughters[], Int_t NDaughters,
681                                    const AliKFParticleBase *Parent,  Double_t Mass         )
682
683   //* Full reconstruction in one go
684
685   Int_t maxIter = 1;
686   bool wasLinearized = fIsLinearized;
687   if( !fIsLinearized ){
688     //fVtxGuess[0] = fVtxGuess[1] = fVtxGuess[2] = 0;  //!!!!
689     fVtxGuess[0] = GetX();
690     fVtxGuess[1] = GetY();
691     fVtxGuess[2] = GetZ();
692     fIsLinearized = 1;
693     maxIter = 3;
694   }
695
696   for( Int_t iter=0; iter<maxIter; iter++ ){
697     fAtProductionVertex = 0;
698     fSFromDecay = 0;
699     fP[0] = fVtxGuess[0];
700     fP[1] = fVtxGuess[1];
701     fP[2] = fVtxGuess[2];
702     fP[3] = 0;
703     fP[4] = 0;
704     fP[5] = 0;
705     fP[6] = 0;
706     fP[7] = 0;
707   
708     for(Int_t i=0;i<36;++i) fC[i]=0.;
709
710     fC[0] = fC[2] = fC[5] = 100.;
711     fC[35] = 1.;
712     
713     fNDF  = -3;
714     fChi2 =  0.;
715     fQ = 0;
716
717     for( Int_t itr =0; itr<NDaughters; itr++ ){
718       AddDaughter( *vDaughters[itr] );    
719     }
720     if( iter<maxIter-1){
721       for( Int_t i=0; i<3; i++ ) fVtxGuess[i] = fP[i];  
722     }
723   }
724   fIsLinearized = wasLinearized;    
725
726   if( Mass>=0 ) SetMassConstraint( Mass );
727   if( Parent ) SetProductionVertex( *Parent );
728 }
729
730
731 void AliKFParticleBase::Convert( bool ToProduction )
732 {
733   //* Tricky function - convert the particle error along its trajectory to 
734   //* the value which corresponds to its production/decay vertex
735   //* It is done by combination of the error of decay length with the position errors
736
737   Double_t fld[3]; 
738   {
739     GetFieldValue( fP, fld );
740     const Double_t kCLight =  fQ*0.000299792458;
741     fld[0]*=kCLight; fld[1]*=kCLight; fld[2]*=kCLight;
742   }
743
744   Double_t h[6];
745   
746   h[0] = fP[3];
747   h[1] = fP[4];
748   h[2] = fP[5];
749   if( ToProduction ){ h[0]=-h[0]; h[1]=-h[1]; h[2]=-h[2]; } 
750   h[3] = h[1]*fld[2]-h[2]*fld[1];
751   h[4] = h[2]*fld[0]-h[0]*fld[2];
752   h[5] = h[0]*fld[1]-h[1]*fld[0];
753   
754   Double_t c;
755
756   c = fC[28]+h[0]*fC[35];
757   fC[ 0]+= h[0]*(c+fC[28]);
758   fC[28] = c;
759
760   fC[ 1]+= h[1]*fC[28] + h[0]*fC[29];
761   c = fC[29]+h[1]*fC[35];
762   fC[ 2]+= h[1]*(c+fC[29]);
763   fC[29] = c;
764
765   fC[ 3]+= h[2]*fC[28] + h[0]*fC[30];
766   fC[ 4]+= h[2]*fC[29] + h[1]*fC[30];
767   c = fC[30]+h[2]*fC[35];
768   fC[ 5]+= h[2]*(c+fC[30]);
769   fC[30] = c;
770
771   fC[ 6]+= h[3]*fC[28] + h[0]*fC[31];
772   fC[ 7]+= h[3]*fC[29] + h[1]*fC[31];
773   fC[ 8]+= h[3]*fC[30] + h[2]*fC[31];
774   c = fC[31]+h[3]*fC[35];
775   fC[ 9]+= h[3]*(c+fC[31]);
776   fC[31] = c;
777   
778   fC[10]+= h[4]*fC[28] + h[0]*fC[32];
779   fC[11]+= h[4]*fC[29] + h[1]*fC[32];
780   fC[12]+= h[4]*fC[30] + h[2]*fC[32];
781   fC[13]+= h[4]*fC[31] + h[3]*fC[32];
782   c = fC[32]+h[4]*fC[35];
783   fC[14]+= h[4]*(c+fC[32]);
784   fC[32] = c;
785   
786   fC[15]+= h[5]*fC[28] + h[0]*fC[33];
787   fC[16]+= h[5]*fC[29] + h[1]*fC[33];
788   fC[17]+= h[5]*fC[30] + h[2]*fC[33];
789   fC[18]+= h[5]*fC[31] + h[3]*fC[33];
790   fC[19]+= h[5]*fC[32] + h[4]*fC[33];
791   c = fC[33]+h[5]*fC[35];
792   fC[20]+= h[5]*(c+fC[33]);
793   fC[33] = c;
794
795   fC[21]+= h[0]*fC[34];
796   fC[22]+= h[1]*fC[34];
797   fC[23]+= h[2]*fC[34];
798   fC[24]+= h[3]*fC[34];
799   fC[25]+= h[4]*fC[34];
800   fC[26]+= h[5]*fC[34];
801 }
802
803
804 void AliKFParticleBase::TransportToDecayVertex()
805 {
806   //* Transport the particle to its decay vertex 
807
808   if( fSFromDecay != 0 ) TransportToDS( -fSFromDecay );
809   if( fAtProductionVertex ) Convert(0);
810   fAtProductionVertex = 0;
811 }
812
813 void AliKFParticleBase::TransportToProductionVertex()
814 {
815   //* Transport the particle to its production vertex 
816   
817   if( fSFromDecay != -fP[7] ) TransportToDS( -fSFromDecay-fP[7] );
818   if( !fAtProductionVertex ) Convert( 1 );
819   fAtProductionVertex = 1;
820 }
821
822
823 void AliKFParticleBase::TransportToDS( Double_t dS )
824
825   //* Transport the particle on dS parameter (SignedPath/Momentum) 
826  
827   Transport( dS, fP, fC );
828   fSFromDecay+= dS;
829 }
830
831
832 Double_t AliKFParticleBase::GetDStoPointLine( const Double_t xyz[] ) const 
833 {
834   //* Get dS to a certain space point without field
835
836   Double_t p2 = fP[3]*fP[3] + fP[4]*fP[4] + fP[5]*fP[5];  
837   if( p2<1.e-4 ) p2 = 1;
838   return ( fP[3]*(xyz[0]-fP[0]) + fP[4]*(xyz[1]-fP[1]) + fP[5]*(xyz[2]-fP[2]) )/p2;
839 }
840
841
842 Double_t AliKFParticleBase::GetDStoPointBz( Double_t B, const Double_t xyz[] ) 
843   const
844
845   
846   //* Get dS to a certain space point for Bz field
847   const Double_t kCLight = 0.000299792458;
848   Double_t bq = B*fQ*kCLight;
849   Double_t pt2 = fP[3]*fP[3] + fP[4]*fP[4];
850   if( pt2<1.e-4 ) return 0;
851   Double_t dx = xyz[0] - fP[0];
852   Double_t dy = xyz[1] - fP[1]; 
853   Double_t a = dx*fP[3]+dy*fP[4];
854   if( TMath::Abs(bq)<1.e-8 ) return a/pt2;  
855   return  TMath::ATan2( bq*a, pt2 + bq*(dy*fP[3] -dx*fP[4]) )/bq;
856 }
857
858
859 void AliKFParticleBase::GetDStoParticleBz( Double_t B, const AliKFParticleBase &p, 
860                                            Double_t &DS, Double_t &DS1 ) 
861   const
862
863   //* Get dS to another particle for Bz field
864   Double_t px = fP[3];
865   Double_t py = fP[4];
866   Double_t pz = fP[5];
867
868   Double_t px1 = p.fP[3];
869   Double_t py1 = p.fP[4];
870   Double_t pz1 = p.fP[5];
871
872   const Double_t kCLight = 0.000299792458;
873
874   Double_t bq = B*fQ*kCLight;
875   Double_t bq1 = B*p.fQ*kCLight;
876   Double_t s=0, ds=0, s1=0, ds1=0;
877   
878   if( TMath::Abs(bq)>1.e-8 || TMath::Abs(bq1)>1.e-8 ){
879
880     Double_t dx = (p.fP[0] - fP[0]);
881     Double_t dy = (p.fP[1] - fP[1]);
882     Double_t d2 = (dx*dx+dy*dy);
883     
884     Double_t p2  = (px *px  + py *py); 
885     Double_t p21 = (px1*px1 + py1*py1);
886     
887     Double_t a = (px*py1 - py*px1);
888     Double_t b = (px*px1 + py*py1);
889     
890     Double_t ldx = bq*bq1*dx - bq1*py + bq*py1 ;
891     Double_t ldy = bq*bq1*dy + bq1*px - bq*px1 ;
892     Double_t l2 = ldx*ldx + ldy*ldy;
893     
894     Double_t cS = bq1*p2 + bq*bq1*(dy* px - dx* py) -  bq*b;
895     Double_t cS1= bq*p21 - bq*bq1*(dy*px1 - dx*py1) - bq1*b;
896
897     Double_t ca  = bq*bq*bq1*d2  +2*( cS + bq*bq*(py1*dx-px1*dy)) ;
898     Double_t ca1 = bq*bq1*bq1*d2 +2*( cS1 - bq1*bq1*(py*dx-px*dy)) ;  
899   
900     Double_t sa = 4*l2*p2 - ca*ca;
901     Double_t sa1 = 4*l2*p21 - ca1*ca1;
902
903     if(sa<0) sa=0;
904     if(sa1<0)sa1=0;
905
906     if( TMath::Abs(bq)>1.e-8){
907       s  = TMath::ATan2(   bq*( bq1*(dx*px +dy*py) + a ) , cS )/bq;
908       ds = TMath::ATan2(TMath::Sqrt(sa),ca)/bq;
909     } else {
910       s = ( (dx*px + dy*py) + (py*px1-px*py1)/bq1)/p2;
911       ds = s*s - (d2-2*(px1*dy-py1*dx)/bq1)/p2; 
912       if( ds<0 ) ds = 0;
913       ds = TMath::Sqrt(ds);   
914     }
915     
916     if( TMath::Abs(bq1)>1.e-8){
917       s1 = TMath::ATan2( -bq1*( bq*(dx*px1+dy*py1) + a), cS1 )/bq1;
918       ds1 = TMath::ATan2(TMath::Sqrt(sa1),ca1)/bq1;  
919     } else {
920       s1 = (-(dx*px1 + dy*py1) + (py*px1-px*py1)/bq)/p21;
921       ds1 = s1*s1 - (d2+2*(px*dy-py*dx)/bq)/p21; 
922       if( ds1<0 ) ds1 = 0;
923       ds1 = TMath::Sqrt(ds1);
924     }
925   }
926
927   Double_t ss[2], ss1[2], g[2][5],g1[2][5];
928   
929   ss[0] = s + ds;
930   ss[1] = s - ds;
931   ss1[0] = s1 + ds1;
932   ss1[1] = s1 - ds1;
933   for( Int_t i=0; i<2; i++){
934     Double_t bs = bq*ss[i];
935     Double_t c = TMath::Cos(bs), s3 = TMath::Sin(bs);
936     Double_t cB,sB;
937     if( TMath::Abs(bq)>1.e-8){
938       cB= (1-c)/bq;     
939       sB= s3/bq;  
940     }else{
941       sB = (1. - bs*bs/6.)*ss[i];
942       cB = .5*sB*bs;
943     }
944     g[i][0] = fP[0] + sB*px + cB*py;
945     g[i][1] = fP[1] - cB*px + sB*py;
946     g[i][2] = fP[2] + ss[i]*pz;
947     g[i][3] =       + c*px + s3*py;
948     g[i][4] =       - s3*px + c*py;
949
950     bs = bq1*ss1[i];  
951     c =  TMath::Cos(bs); s3 = TMath::Sin(bs);
952     if( TMath::Abs(bq1)>1.e-8){
953       cB= (1-c)/bq1;   
954       sB= s3/bq1;  
955     }else{
956       sB = (1. - bs*bs/6.)*ss1[i];
957       cB = .5*sB*bs;
958     }
959       
960     g1[i][0] = p.fP[0] + sB*px1 + cB*py1;
961     g1[i][1] = p.fP[1] - cB*px1 + sB*py1;
962     g1[i][2] = p.fP[2] + ss[i]*pz1;
963     g1[i][3] =         + c*px1 + s3*py1;
964     g1[i][4] =         - s3*px1 + c*py1;
965   }
966
967   Int_t i=0, i1=0;
968   
969   Double_t dMin = 1.e10;
970   for( Int_t j=0; j<2; j++){
971     for( Int_t j1=0; j1<2; j1++){
972       Double_t xx = g[j][0]-g1[j1][0];
973       Double_t yy = g[j][1]-g1[j1][1];
974       Double_t zz = g[j][2]-g1[j1][2];
975       Double_t d = xx*xx + yy*yy + zz*zz;
976       if( d<dMin ){
977         dMin = d;
978         i = j;
979         i1 = j1;
980       }
981     }
982   }  
983
984   DS = ss[i];
985   DS1 = ss1[i1];
986   if(0){
987     Double_t x= g[i][0], y= g[i][1], z= g[i][2], ppx= g[i][3], ppy= g[i][4];  
988     Double_t x1=g1[i1][0], y1= g1[i1][1], z1= g1[i1][2], ppx1= g1[i1][3], ppy1= g1[i1][4];  
989     Double_t dx = x1-x;
990     Double_t dy = y1-y;
991     Double_t dz = z1-z;
992     Double_t a = ppx*ppx1 + ppy*ppy1 + pz*pz1;
993     Double_t b = dx*ppx1 + dy*ppy1 + dz*pz1;
994     Double_t c = dx*ppx  + dy*ppy  + dz*pz ;
995     Double_t pp2 = ppx*ppx + ppy*ppy + pz*pz;
996     Double_t pp21= ppx1*ppx1 + ppy1*ppy1 + pz1*pz1;
997     Double_t det = pp2*pp21 - a*a;
998     if( TMath::Abs(det)>1.e-8 ){
999       DS+=(a*b-pp21*c)/det;
1000       DS1+=(a*c-pp2*b)/det;
1001     }
1002   }
1003 }
1004
1005
1006
1007 void AliKFParticleBase::TransportCBM( Double_t dS, 
1008                                  Double_t P[], Double_t C[] ) const
1009 {  
1010   //* Transport the particle on dS, output to P[],C[], for CBM field
1011  
1012   if( fQ==0 ){
1013     TransportLine( dS, P, C );
1014     return;
1015   }
1016
1017   const Double_t kCLight = 0.000299792458;
1018
1019   Double_t c = fQ*kCLight;
1020
1021   // construct coefficients 
1022
1023   Double_t 
1024     px   = fP[3],
1025     py   = fP[4],
1026     pz   = fP[5];
1027       
1028   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;
1029
1030   { // get field integrals
1031
1032     Double_t fld[3][3];   
1033     Double_t p0[3], p1[3], p2[3];
1034
1035     // line track approximation
1036
1037     p0[0] = fP[0];
1038     p0[1] = fP[1];
1039     p0[2] = fP[2];
1040   
1041     p2[0] = fP[0] + px*dS;
1042     p2[1] = fP[1] + py*dS;
1043     p2[2] = fP[2] + pz*dS;
1044   
1045     p1[0] = 0.5*(p0[0]+p2[0]);
1046     p1[1] = 0.5*(p0[1]+p2[1]);
1047     p1[2] = 0.5*(p0[2]+p2[2]);
1048
1049     // first order track approximation
1050     {
1051       GetFieldValue( p0, fld[0] );
1052       GetFieldValue( p1, fld[1] );
1053       GetFieldValue( p2, fld[2] );
1054
1055       Double_t ssy1 = ( 7*fld[0][1] + 6*fld[1][1]-fld[2][1] )*c*dS*dS/96.;
1056       Double_t ssy2 = (   fld[0][1] + 2*fld[1][1]         )*c*dS*dS/6.;
1057
1058       p1[0] -= ssy1*pz;
1059       p1[2] += ssy1*px;
1060       p2[0] -= ssy2*pz;
1061       p2[2] += ssy2*px;   
1062     }
1063
1064     GetFieldValue( p0, fld[0] );
1065     GetFieldValue( p1, fld[1] );
1066     GetFieldValue( p2, fld[2] );
1067     
1068     sx = c*( fld[0][0] + 4*fld[1][0] + fld[2][0] )*dS/6.;
1069     sy = c*( fld[0][1] + 4*fld[1][1] + fld[2][1] )*dS/6.;
1070     sz = c*( fld[0][2] + 4*fld[1][2] + fld[2][2] )*dS/6.;
1071
1072     ssx = c*( fld[0][0] + 2*fld[1][0])*dS*dS/6.;
1073     ssy = c*( fld[0][1] + 2*fld[1][1])*dS*dS/6.;
1074     ssz = c*( fld[0][2] + 2*fld[1][2])*dS*dS/6.;
1075
1076     Double_t c2[3][3]    =   { {  5, -4, -1},{  44,  80,  -4},{ 11, 44, 5} }; // /=360.    
1077     Double_t cc2[3][3]    =   { { 38,  8, -4},{ 148, 208, -20},{  3, 36, 3} }; // /=2520.
1078     for(Int_t n=0; n<3; n++)
1079       for(Int_t m=0; m<3; m++) 
1080         {
1081           syz += c2[n][m]*fld[n][1]*fld[m][2];
1082           ssyz += cc2[n][m]*fld[n][1]*fld[m][2];
1083         }
1084  
1085     syz  *= c*c*dS*dS/360.;
1086     ssyz  *= c*c*dS*dS*dS/2520.;
1087     
1088     syy  = c*( fld[0][1] + 4*fld[1][1] + fld[2][1] )*dS;
1089     syyy = syy*syy*syy / 1296;
1090     syy  = syy*syy/72;
1091
1092     ssyy = ( fld[0][1]*( 38*fld[0][1] + 156*fld[1][1]  -   fld[2][1] )+
1093             fld[1][1]*(              208*fld[1][1]  +16*fld[2][1] )+
1094             fld[2][1]*(                             3*fld[2][1] )  
1095             )*dS*dS*dS*c*c/2520.;
1096     ssyyy = 
1097       (
1098        fld[0][1]*( fld[0][1]*( 85*fld[0][1] + 526*fld[1][1]  - 7*fld[2][1] )+
1099                  fld[1][1]*(             1376*fld[1][1]  +84*fld[2][1] )+
1100                  fld[2][1]*(                            19*fld[2][1] )  )+
1101        fld[1][1]*( fld[1][1]*(             1376*fld[1][1] +256*fld[2][1] )+
1102                  fld[2][1]*(                            62*fld[2][1] )  )+
1103        fld[2][1]*fld[2][1]  *(                             3*fld[2][1] )       
1104        )*dS*dS*dS*dS*c*c*c/90720.;    
1105  
1106   }
1107
1108   Double_t mJ[8][8];
1109   for( Int_t i=0; i<8; i++ ) for( Int_t j=0; j<8; j++) mJ[i][j]=0;
1110
1111   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;
1112   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;
1113   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;
1114   
1115   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;
1116   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;
1117   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;
1118   mJ[6][6] = mJ[7][7] = 1;
1119   
1120   P[0] = fP[0] + mJ[0][3]*px + mJ[0][4]*py + mJ[0][5]*pz;
1121   P[1] = fP[1] + mJ[1][3]*px + mJ[1][4]*py + mJ[1][5]*pz;
1122   P[2] = fP[2] + mJ[2][3]*px + mJ[2][4]*py + mJ[2][5]*pz;
1123   P[3] =        mJ[3][3]*px + mJ[3][4]*py + mJ[3][5]*pz;
1124   P[4] =        mJ[4][3]*px + mJ[4][4]*py + mJ[4][5]*pz;
1125   P[5] =        mJ[5][3]*px + mJ[5][4]*py + mJ[5][5]*pz;
1126   P[6] = fP[6];
1127   P[7] = fP[7];
1128
1129   MultQSQt( mJ[0], fC, C);
1130
1131 }
1132
1133
1134 void AliKFParticleBase::TransportBz( Double_t B, Double_t dS,
1135                                      Double_t P[], Double_t C[] ) const 
1136
1137   //* Transport the particle on dS, output to P[],C[], for Bz field
1138  
1139   const Double_t kCLight = 0.000299792458;
1140   B = B*fQ*kCLight;
1141   Double_t bs= B*dS;
1142   Double_t s = TMath::Sin(bs), c = TMath::Cos(bs);
1143   Double_t sB, cB;
1144   if( TMath::Abs(bs)>1.e-10){
1145     sB= s/B;
1146     cB= (1-c)/B;
1147   }else{
1148     sB = (1. - bs*bs/6.)*dS;
1149     cB = .5*sB*bs;
1150   }
1151   
1152   Double_t px = fP[3];
1153   Double_t py = fP[4];
1154   Double_t pz = fP[5];
1155   
1156   P[0] = fP[0] + sB*px + cB*py;
1157   P[1] = fP[1] - cB*px + sB*py;
1158   P[2] = fP[2] +  dS*pz;
1159   P[3] =          c*px + s*py;
1160   P[4] =         -s*px + c*py;
1161   P[5] = fP[5];
1162   P[6] = fP[6];
1163   P[7] = fP[7];
1164  
1165   Double_t mJ[8][8] = { {1,0,0,   sB, cB,  0, 0, 0 },
1166                         {0,1,0,  -cB, sB,  0, 0, 0 },
1167                         {0,0,1,    0,  0,  dS, 0, 0 },
1168                         {0,0,0,    c,  s,  0, 0, 0 },
1169                         {0,0,0,   -s,  c,  0, 0, 0 },
1170                         {0,0,0,    0,  0,  1, 0, 0 },
1171                         {0,0,0,    0,  0,  0, 1, 0 },
1172                         {0,0,0,    0,  0,  0, 0, 1 }  };
1173   Double_t mA[8][8];
1174   for( Int_t k=0,i=0; i<8; i++)
1175     for( Int_t j=0; j<=i; j++, k++ ) mA[i][j] = mA[j][i] = fC[k]; 
1176
1177   Double_t mJC[8][8];
1178   for( Int_t i=0; i<8; i++ )
1179     for( Int_t j=0; j<8; j++ ){
1180       mJC[i][j]=0;
1181       for( Int_t k=0; k<8; k++ ) mJC[i][j]+=mJ[i][k]*mA[k][j];
1182     }
1183   
1184   for( Int_t k=0,i=0; i<8; i++)
1185     for( Int_t j=0; j<=i; j++, k++ ){
1186       C[k] = 0;
1187       for( Int_t l=0; l<8; l++ ) C[k]+=mJC[i][l]*mJ[j][l];
1188     }
1189   
1190   return;
1191   /*
1192   Double_t cBC13 = cB*fC[13];
1193   Double_t C17 = fC[17];
1194   Double_t C18 = fC[18];
1195   fC[ 0]+=  2*(sB*fC[ 6] + cB*fC[10]) + (sB*fC[ 9] + 2*cBC13)*sB + cB*cB*fC[14];
1196
1197   Double_t mJC13= fC[ 7] - cB*fC[ 9] + sB*fC[13];
1198   Double_t mJC14= fC[11] -    cBC13 + sB*fC[14];
1199  
1200   fC[ 1]+= -cB*fC[ 6] + sB*fC[10] +mJC13*sB +mJC14*cB;
1201   fC[ 2]+= -cB*fC[ 7] + sB*fC[11] -mJC13*cB +mJC14*sB;
1202
1203   Double_t mJC23= fC[ 8] + dS*fC[18];
1204   Double_t mJC24= fC[12] + dS*fC[19];
1205   fC[ 3]+= dS*fC[15] +mJC23*sB +mJC24*cB;
1206   fC[ 4]+= dS*fC[16] -mJC23*cB +mJC24*sB;
1207  
1208   fC[15]+=  C18*sB + fC[19]*cB;
1209   fC[16]+= -C18*cB + fC[19]*sB;
1210   fC[17]+=  fC[20]*dS;
1211   fC[18] =  C18*c + fC[19]*s;
1212   fC[19] = -C18*s + fC[19]*c;
1213
1214   fC[ 5]+= (C17 + C17 + fC[17])*dS;
1215
1216   Double_t mJC33= c*fC[ 9] + s*fC[13]; Double_t mJC34= c*fC[13] + s*fC[14];
1217   Double_t mJC43=-s*fC[ 9] + c*fC[13]; Double_t mJC44=-s*fC[13] + c*fC[14];
1218   Double_t C6= fC[6], C7= fC[7], C8= fC[8]; 
1219
1220   fC[ 6]= c*C6 + s*fC[10] +mJC33*sB +mJC34*cB;
1221   fC[ 7]= c*C7 + s*fC[11] -mJC33*cB +mJC34*sB;
1222   fC[ 8]= c*C8 + s*fC[12] +fC[18]*dS;
1223   fC[ 9]= mJC33*c +mJC34*s;
1224
1225   fC[10]= -s*C6 + c*fC[10] +mJC43*sB +mJC44*cB;
1226   fC[11]= -s*C7 + c*fC[11] -mJC43*cB +mJC44*sB;
1227   fC[12]= -s*C8 + c*fC[12] +fC[19]*dS;
1228   fC[13]= mJC43*c +mJC44*s;
1229   fC[14]=-mJC43*s +mJC44*c;
1230   */
1231 }
1232
1233
1234 Double_t AliKFParticleBase::GetDistanceFromVertex( const AliKFParticleBase &Vtx ) const
1235 {
1236   //* Calculate distance from vertex [cm]
1237
1238   return GetDistanceFromVertex( Vtx.fP );
1239 }
1240
1241 Double_t AliKFParticleBase::GetDistanceFromVertex( const Double_t vtx[] ) const
1242 {
1243   //* Calculate distance from vertex [cm]
1244
1245   Double_t mP[8], mC[36];  
1246   Transport( GetDStoPoint(vtx), mP, mC );
1247   Double_t d[3]={ vtx[0]-mP[0], vtx[1]-mP[1], vtx[2]-mP[2]};
1248   return TMath::Sqrt( d[0]*d[0]+d[1]*d[1]+d[2]*d[2] );
1249 }
1250
1251 Double_t AliKFParticleBase::GetDistanceFromParticle( const AliKFParticleBase &p ) 
1252   const
1253
1254   //* Calculate distance to other particle [cm]
1255
1256   Double_t dS, dS1;
1257   GetDStoParticle( p, dS, dS1 );   
1258   Double_t mP[8], mC[36], mP1[8], mC1[36];
1259   Transport( dS, mP, mC ); 
1260   p.Transport( dS1, mP1, mC1 ); 
1261   Double_t dx = mP[0]-mP1[0]; 
1262   Double_t dy = mP[1]-mP1[1]; 
1263   Double_t dz = mP[2]-mP1[2]; 
1264   return TMath::Sqrt(dx*dx+dy*dy+dz*dz);
1265 }
1266
1267 Double_t AliKFParticleBase::GetDeviationFromVertex( const AliKFParticleBase &Vtx ) const
1268 {
1269   //* Calculate sqrt(Chi2/ndf) deviation from vertex
1270
1271   return GetDeviationFromVertex( Vtx.fP, Vtx.fC );
1272 }
1273
1274
1275 Double_t AliKFParticleBase::GetDeviationFromVertex( const Double_t v[], const Double_t Cv[] ) const
1276 {
1277   //* Calculate sqrt(Chi2/ndf) deviation from vertex
1278   //* v = [xyz], Cv=[Cxx,Cxy,Cyy,Cxz,Cyz,Czz]-covariance matrix
1279
1280   Double_t mP[8];
1281   Double_t mC[36];
1282   
1283   Transport( GetDStoPoint(v), mP, mC );  
1284
1285   Double_t d[3]={ v[0]-mP[0], v[1]-mP[1], v[2]-mP[2]};
1286
1287   Double_t sigmaS = .1+10.*TMath::Sqrt( (d[0]*d[0]+d[1]*d[1]+d[2]*d[2])/
1288                                  (mP[3]*mP[3]+mP[4]*mP[4]+mP[5]*mP[5])  );
1289
1290    
1291   Double_t h[3] = { mP[3]*sigmaS, mP[4]*sigmaS, mP[5]*sigmaS };       
1292   
1293   Double_t mSi[6] = 
1294     { mC[0] +h[0]*h[0], 
1295       mC[1] +h[1]*h[0], mC[2] +h[1]*h[1], 
1296       mC[3] +h[2]*h[0], mC[4] +h[2]*h[1], mC[5] +h[2]*h[2] };
1297
1298   if( Cv ){
1299     mSi[0]+=Cv[0];
1300     mSi[1]+=Cv[1];
1301     mSi[2]+=Cv[2];
1302     mSi[3]+=Cv[3];
1303     mSi[4]+=Cv[4];
1304     mSi[5]+=Cv[5];
1305   }
1306   
1307   Double_t mS[6];
1308
1309   mS[0] = mSi[2]*mSi[5] - mSi[4]*mSi[4];
1310   mS[1] = mSi[3]*mSi[4] - mSi[1]*mSi[5];
1311   mS[2] = mSi[0]*mSi[5] - mSi[3]*mSi[3];
1312   mS[3] = mSi[1]*mSi[4] - mSi[2]*mSi[3];
1313   mS[4] = mSi[1]*mSi[3] - mSi[0]*mSi[4];
1314   mS[5] = mSi[0]*mSi[2] - mSi[1]*mSi[1];         
1315       
1316   Double_t s = ( mSi[0]*mS[0] + mSi[1]*mS[1] + mSi[3]*mS[3] );
1317   s = ( s > 1.E-20 )  ?1./s :0;   
1318
1319   return TMath::Sqrt( TMath::Abs(s*( ( mS[0]*d[0] + mS[1]*d[1] + mS[3]*d[2])*d[0]
1320                    +(mS[1]*d[0] + mS[2]*d[1] + mS[4]*d[2])*d[1]
1321                    +(mS[3]*d[0] + mS[4]*d[1] + mS[5]*d[2])*d[2] ))/2);
1322 }
1323
1324
1325 Double_t AliKFParticleBase::GetDeviationFromParticle( const AliKFParticleBase &p ) 
1326   const
1327
1328   //* Calculate sqrt(Chi2/ndf) deviation from other particle
1329
1330   Double_t dS, dS1;
1331   GetDStoParticle( p, dS, dS1 );   
1332   Double_t mP1[8], mC1[36];
1333   p.Transport( dS1, mP1, mC1 ); 
1334
1335   Double_t d[3]={ fP[0]-mP1[0], fP[1]-mP1[1], fP[2]-mP1[2]};
1336
1337   Double_t sigmaS = .1+10.*TMath::Sqrt( (d[0]*d[0]+d[1]*d[1]+d[2]*d[2])/
1338                                         (mP1[3]*mP1[3]+mP1[4]*mP1[4]+mP1[5]*mP1[5])  );
1339
1340   Double_t h[3] = { mP1[3]*sigmaS, mP1[4]*sigmaS, mP1[5]*sigmaS };       
1341   
1342   mC1[0] +=h[0]*h[0];
1343   mC1[1] +=h[1]*h[0]; 
1344   mC1[2] +=h[1]*h[1]; 
1345   mC1[3] +=h[2]*h[0]; 
1346   mC1[4] +=h[2]*h[1];
1347   mC1[5] +=h[2]*h[2];
1348
1349   return GetDeviationFromVertex( mP1, mC1 )*TMath::Sqrt(2./1.);
1350 }
1351
1352
1353
1354 void AliKFParticleBase::SubtractFromVertex( Double_t v[], Double_t Cv[], 
1355                                      Double_t &vChi2, Int_t vNDF ) const 
1356 {    
1357   //* Subtract the particle from the vertex  
1358   
1359   Double_t fld[3];  
1360   {
1361     GetFieldValue( v, fld );
1362     const Double_t kCLight =  0.000299792458;
1363     fld[0]*=kCLight; fld[1]*=kCLight; fld[2]*=kCLight;
1364   }
1365
1366   Double_t m[8];
1367   Double_t mCm[36];
1368     
1369   Transport( GetDStoPoint(v), m, mCm );
1370
1371   Double_t d[3] = { v[0]-m[0], v[1]-m[1], v[2]-m[2] };
1372   Double_t sigmaS = .1+10.*TMath::Sqrt( (d[0]*d[0]+d[1]*d[1]+d[2]*d[2])/
1373                                         (m[3]*m[3]+m[4]*m[4]+m[5]*m[5])  );
1374   
1375   Double_t h[6];
1376
1377   h[0] = m[3]*sigmaS;
1378   h[1] = m[4]*sigmaS;
1379   h[2] = m[5]*sigmaS; 
1380   h[3] = ( h[1]*fld[2]-h[2]*fld[1] )*GetQ();
1381   h[4] = ( h[2]*fld[0]-h[0]*fld[2] )*GetQ();
1382   h[5] = ( h[0]*fld[1]-h[1]*fld[0] )*GetQ();
1383     
1384   //* Fit of daughter momentum (Px,Py,Pz) to fVtxGuess vertex
1385   {
1386     Double_t zeta[3] = { v[0]-m[0], v[1]-m[1], v[2]-m[2] };
1387     
1388     Double_t mVv[6] = 
1389       { mCm[ 0] + h[0]*h[0], 
1390         mCm[ 1] + h[1]*h[0], mCm[ 2] + h[1]*h[1],                            
1391         mCm[ 3] + h[2]*h[0], mCm[ 4] + h[2]*h[1], mCm[ 5] + h[2]*h[2] };
1392     
1393     Double_t mVvp[9]=
1394       { mCm[ 6] + h[0]*h[3], mCm[ 7] + h[1]*h[3], mCm[ 8] + h[2]*h[3],
1395         mCm[10] + h[0]*h[4], mCm[11] + h[1]*h[4], mCm[12] + h[2]*h[4],
1396         mCm[15] + h[0]*h[5], mCm[16] + h[1]*h[5], mCm[17] + h[2]*h[5] };
1397       
1398     Double_t mS[6] = 
1399       { mVv[2]*mVv[5] - mVv[4]*mVv[4],
1400         mVv[3]*mVv[4] - mVv[1]*mVv[5], mVv[0]*mVv[5] - mVv[3]*mVv[3],
1401         mVv[1]*mVv[4] - mVv[2]*mVv[3], mVv[1]*mVv[3] - mVv[0]*mVv[4], 
1402         mVv[0]*mVv[2] - mVv[1]*mVv[1] };                 
1403       
1404     Double_t s = ( mVv[0]*mS[0] + mVv[1]*mS[1] + mVv[3]*mS[3] );
1405     s = ( s > 1.E-20 )  ?1./s :0;
1406     
1407     mS[0]*=s; mS[1]*=s; mS[2]*=s; mS[3]*=s; mS[4]*=s; mS[5]*=s;
1408       
1409     Double_t mSz[3] = { (mS[0]*zeta[0]+mS[1]*zeta[1]+mS[3]*zeta[2]),
1410                        (mS[1]*zeta[0]+mS[2]*zeta[1]+mS[4]*zeta[2]),
1411                        (mS[3]*zeta[0]+mS[4]*zeta[1]+mS[5]*zeta[2]) };
1412     
1413     Double_t px = m[3] + mVvp[0]*mSz[0] + mVvp[1]*mSz[1] + mVvp[2]*mSz[2];
1414     Double_t py = m[4] + mVvp[3]*mSz[0] + mVvp[4]*mSz[1] + mVvp[5]*mSz[2];
1415     Double_t pz = m[5] + mVvp[6]*mSz[0] + mVvp[7]*mSz[1] + mVvp[8]*mSz[2];
1416     
1417     h[0] = px*sigmaS;
1418     h[1] = py*sigmaS;
1419     h[2] = pz*sigmaS; 
1420     h[3] = ( h[1]*fld[2]-h[2]*fld[1] )*GetQ();
1421     h[4] = ( h[2]*fld[0]-h[0]*fld[2] )*GetQ();
1422     h[5] = ( h[0]*fld[1]-h[1]*fld[0] )*GetQ();
1423   }
1424     
1425   Double_t mV[6];
1426     
1427   mV[ 0] = mCm[ 0] + h[0]*h[0];
1428   mV[ 1] = mCm[ 1] + h[1]*h[0];
1429   mV[ 2] = mCm[ 2] + h[1]*h[1];
1430   mV[ 3] = mCm[ 3] + h[2]*h[0];
1431   mV[ 4] = mCm[ 4] + h[2]*h[1];
1432   mV[ 5] = mCm[ 5] + h[2]*h[2];
1433      
1434   //* 
1435             
1436   Double_t mS[6];
1437   {
1438     Double_t mSi[6] = { mV[0]-Cv[0], 
1439                        mV[1]-Cv[1], mV[2]-Cv[2], 
1440                        mV[3]-Cv[3], mV[4]-Cv[4], mV[5]-Cv[5] };
1441     
1442     mS[0] = mSi[2]*mSi[5] - mSi[4]*mSi[4];
1443     mS[1] = mSi[3]*mSi[4] - mSi[1]*mSi[5];
1444     mS[2] = mSi[0]*mSi[5] - mSi[3]*mSi[3];
1445     mS[3] = mSi[1]*mSi[4] - mSi[2]*mSi[3];
1446     mS[4] = mSi[1]*mSi[3] - mSi[0]*mSi[4];
1447     mS[5] = mSi[0]*mSi[2] - mSi[1]*mSi[1];       
1448     
1449     Double_t s = ( mSi[0]*mS[0] + mSi[1]*mS[1] + mSi[3]*mS[3] );
1450     s = ( s > 1.E-20 )  ?1./s :0;         
1451     mS[0]*=s;
1452     mS[1]*=s;
1453     mS[2]*=s;
1454     mS[3]*=s;
1455     mS[4]*=s;
1456     mS[5]*=s;
1457   }
1458     
1459   //* Residual (measured - estimated)
1460     
1461   Double_t zeta[3] = { m[0]-v[0], m[1]-v[1], m[2]-v[2] };
1462         
1463   //* mCHt = mCH' - D'
1464     
1465   Double_t mCHt0[3], mCHt1[3], mCHt2[3];
1466     
1467   mCHt0[0]=Cv[ 0] ;      mCHt1[0]=Cv[ 1] ;      mCHt2[0]=Cv[ 3] ;
1468   mCHt0[1]=Cv[ 1] ;      mCHt1[1]=Cv[ 2] ;      mCHt2[1]=Cv[ 4] ;
1469   mCHt0[2]=Cv[ 3] ;      mCHt1[2]=Cv[ 4] ;      mCHt2[2]=Cv[ 5] ;
1470   
1471   //* Kalman gain K = mCH'*S
1472     
1473   Double_t k0[3], k1[3], k2[3];
1474     
1475   for(Int_t i=0;i<3;++i){
1476     k0[i] = mCHt0[i]*mS[0] + mCHt1[i]*mS[1] + mCHt2[i]*mS[3];
1477     k1[i] = mCHt0[i]*mS[1] + mCHt1[i]*mS[2] + mCHt2[i]*mS[4];
1478     k2[i] = mCHt0[i]*mS[3] + mCHt1[i]*mS[4] + mCHt2[i]*mS[5];
1479   }
1480     
1481   //* New estimation of the vertex position r += K*zeta
1482     
1483   for(Int_t i=0;i<3;++i) 
1484     v[i] -= k0[i]*zeta[0] + k1[i]*zeta[1] + k2[i]*zeta[2];
1485     
1486   //* New covariance matrix C -= K*(mCH')'
1487     
1488   for(Int_t i=0, k=0;i<3;++i){
1489     for(Int_t j=0;j<=i;++j,++k) 
1490       Cv[k] += k0[i]*mCHt0[j] + k1[i]*mCHt1[j] + k2[i]*mCHt2[j];
1491   }
1492     
1493   //* Calculate Chi^2 
1494   
1495   vNDF  -= 2;
1496   vChi2 -= (mS[0]*zeta[0] + mS[1]*zeta[1] + mS[3]*zeta[2])*zeta[0]
1497     +      (mS[1]*zeta[0] + mS[2]*zeta[1] + mS[4]*zeta[2])*zeta[1]
1498     +      (mS[3]*zeta[0] + mS[4]*zeta[1] + mS[5]*zeta[2])*zeta[2];  
1499 }
1500
1501
1502
1503 void AliKFParticleBase::TransportLine( Double_t dS, 
1504                                        Double_t P[], Double_t C[] ) const 
1505 {
1506   //* Transport the particle as a straight line
1507
1508   P[0] = fP[0] + dS*fP[3];
1509   P[1] = fP[1] + dS*fP[4];
1510   P[2] = fP[2] + dS*fP[5];
1511   P[3] = fP[3];
1512   P[4] = fP[4];
1513   P[5] = fP[5];
1514   P[6] = fP[6];
1515   P[7] = fP[7];
1516  
1517   Double_t c6  = fC[ 6] + dS*fC[ 9];
1518   Double_t c11 = fC[11] + dS*fC[14];
1519   Double_t c17 = fC[17] + dS*fC[20];
1520   Double_t sc13 = dS*fC[13];
1521   Double_t sc18 = dS*fC[18];
1522   Double_t sc19 = dS*fC[19];
1523
1524   C[ 0] = fC[ 0] + dS*( fC[ 6] + c6  );
1525   C[ 2] = fC[ 2] + dS*( fC[11] + c11 );
1526   C[ 5] = fC[ 5] + dS*( fC[17] + c17 );
1527
1528   C[ 7] = fC[ 7] + sc13;
1529   C[ 8] = fC[ 8] + sc18;
1530   C[ 9] = fC[ 9];
1531
1532   C[12] = fC[12] + sc19;
1533
1534   C[ 1] = fC[ 1] + dS*( fC[10] + C[ 7] );
1535   C[ 3] = fC[ 3] + dS*( fC[15] + C[ 8] );
1536   C[ 4] = fC[ 4] + dS*( fC[16] + C[12] ); 
1537   C[ 6] = c6;
1538
1539   C[10] = fC[10] + sc13;
1540   C[11] = c11;
1541
1542   C[13] = fC[13];
1543   C[14] = fC[14];
1544   C[15] = fC[15] + sc18;
1545   C[16] = fC[16] + sc19;
1546   C[17] = c17;
1547   
1548   C[18] = fC[18];
1549   C[19] = fC[19];
1550   C[20] = fC[20];
1551   C[21] = fC[21] + dS*fC[24];
1552   C[22] = fC[22] + dS*fC[25];
1553   C[23] = fC[23] + dS*fC[26];
1554
1555   C[24] = fC[24];
1556   C[25] = fC[25];
1557   C[26] = fC[26];
1558   C[27] = fC[27];
1559   C[28] = fC[28] + dS*fC[31];
1560   C[29] = fC[29] + dS*fC[32];
1561   C[30] = fC[30] + dS*fC[33];
1562
1563   C[31] = fC[31];
1564   C[32] = fC[32];
1565   C[33] = fC[33];
1566   C[34] = fC[34];
1567   C[35] = fC[35]; 
1568 }
1569
1570
1571 void AliKFParticleBase::MultQSQt( const Double_t Q[], const Double_t S[], Double_t SOut[] )
1572 {
1573   //* Matrix multiplication Q*S*Q^T, Q - square matrix, S - symmetric
1574
1575   const Int_t kN= 8;
1576   Double_t mA[kN*kN];
1577   
1578   for( Int_t i=0, ij=0; i<kN; i++ ){
1579     for( Int_t j=0; j<kN; j++, ++ij ){
1580       mA[ij] = 0 ;
1581       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];
1582     }
1583   }
1584     
1585   for( Int_t i=0; i<kN; i++ ){
1586     for( Int_t j=0; j<=i; j++ ){
1587       Int_t ij = ( j<=i ) ? i*(i+1)/2+j :j*(j+1)/2+i;
1588       SOut[ij] = 0 ;
1589       for( Int_t k=0; k<kN; k++ )  SOut[ij] += Q[ i*kN+k ] * mA[ k*kN+j ];
1590     }
1591   }
1592 }
1593
1594
1595 // 72-charachters line to define the printer border
1596 //3456789012345678901234567890123456789012345678901234567890123456789012