]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliKFParticleBase.cxx
fixed memory leak
[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 &p, Double_t &error )  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   p = TMath::Sqrt(p2);
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;
121     return 0;
122   }
123   return 1;
124 }
125
126 Int_t AliKFParticleBase::GetPt( Double_t &pt, Double_t &error )  const 
127 {
128   //* Calculate particle transverse momentum
129
130   Double_t px = fP[3];
131   Double_t py = fP[4];
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;
139     return 0;
140   }
141   error = 1.e10;
142   return 1;
143 }
144
145 Int_t AliKFParticleBase::GetEta( Double_t &eta, Double_t &error )  const 
146 {
147   //* Calculate particle pseudorapidity
148
149   Double_t px = fP[3];
150   Double_t py = fP[4];
151   Double_t pz = fP[5];
152   Double_t pt2 = px*px + py*py;
153   Double_t p2 = pt2 + pz*pz;
154   Double_t p = TMath::Sqrt(p2);
155   Double_t a = p + pz;
156   Double_t b = p - pz;
157   eta = 1.e10;
158   if( b > 1.e-8 ){
159     Double_t c = a/b;
160     if( c>1.e-8 ) eta = 0.5*TMath::Log(c);
161   }
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] ) );
167
168   if( error>0 && p2pt4>1.e-10 ){
169     error = TMath::Sqrt(error/p2pt4);
170     return 0;
171   }
172
173   error = 1.e10;
174   return 1;
175 }
176
177 Int_t AliKFParticleBase::GetPhi( Double_t &phi, Double_t &error )  const 
178 {
179   //* Calculate particle polar angle
180
181   Double_t px = fP[3];
182   Double_t py = fP[4];
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;
190     return 0;
191   }
192   error = 1.e10;
193   return 1;
194 }
195
196 Int_t AliKFParticleBase::GetR( Double_t &r, Double_t &error )  const 
197 {
198   //* Calculate distance to the origin
199
200   Double_t x = fP[0];
201   Double_t y = fP[1];
202   Double_t x2 = x*x;
203   Double_t y2 = y*y;
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;
208     return 0;
209   }
210   error = 1.e10;
211   return 1;
212 }
213
214 Int_t AliKFParticleBase::GetMass( Double_t &m, Double_t &error ) const 
215 {
216   //* Calculate particle mass
217   
218   // s = sigma^2 of m2/2
219
220   Double_t s = (  fP[3]*fP[3]*fC[9] + fP[4]*fP[4]*fC[14] + fP[5]*fP[5]*fC[20] 
221                   + fP[6]*fP[6]*fC[27] 
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] )   )
224                  ); 
225   Double_t m2 = TMath::Abs(fP[6]*fP[6] - fP[3]*fP[3] - fP[4]*fP[4] - fP[5]*fP[5]);
226   m  = TMath::Sqrt(m2);
227   if( m>1.e-10 ){
228     if( s>=0 ){
229       error = TMath::Sqrt(s)/m;
230       return 0;
231     }
232   }
233   error = 1.e20;
234   return 1;
235 }
236
237
238 Int_t AliKFParticleBase::GetDecayLength( Double_t &l, Double_t &error ) const 
239 {
240   //* Calculate particle decay length [cm]
241
242   Double_t x = fP[3];
243   Double_t y = fP[4];
244   Double_t z = fP[5];
245   Double_t t = fP[7];
246   Double_t x2 = x*x;
247   Double_t y2 = y*y;
248   Double_t z2 = z*z;
249   Double_t p2 = x2+y2+z2;
250   l = t*TMath::Sqrt(p2);
251   if( p2>1.e-4){
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));
256     return 0;
257   }
258   error = 1.e20;
259   return 1;
260 }
261
262 Int_t AliKFParticleBase::GetDecayLengthXY( Double_t &l, Double_t &error ) const 
263 {
264   //* Calculate particle decay length in XY projection [cm]
265
266   Double_t x = fP[3];
267   Double_t y = fP[4];
268   Double_t t = fP[7];
269   Double_t x2 = x*x;
270   Double_t y2 = y*y;
271   Double_t pt2 = x2+y2;
272   l = t*TMath::Sqrt(pt2);
273   if( pt2>1.e-4){
274     error = pt2*fC[35] + t*t/pt2*(x2*fC[9]+y2*fC[14] + 2*x*y*fC[13] )
275       + 2*t*(x*fC[31]+y*fC[32]);
276     error = TMath::Sqrt(TMath::Abs(error));
277     return 0;
278   }
279   error = 1.e20;
280   return 1;
281 }
282
283
284 Int_t AliKFParticleBase::GetLifeTime( Double_t &tauC, Double_t &error ) const 
285 {
286   //* Calculate particle decay time [s]
287
288   Double_t m, dm;
289   GetMass( m, dm );
290   Double_t cTM = (-fP[3]*fC[31] - fP[4]*fC[32] - fP[5]*fC[33] + fP[6]*fC[34]);
291   tauC = fP[7]*m;
292   error = m*m*fC[35] + 2*fP[7]*cTM + fP[7]*fP[7]*dm*dm;
293   if( error > 0 ){
294     error = TMath::Sqrt( error );
295     return 0;
296   }
297   error = 1.e20;
298   return 1;
299 }
300
301
302 void AliKFParticleBase::operator +=( const AliKFParticleBase &Daughter )
303 {
304   //* Add daughter via operator+=
305
306   AddDaughter( Daughter );
307 }
308   
309 Double_t AliKFParticleBase::GetSCorrection( const Double_t Part[], const Double_t XYZ[] ) 
310 {
311   //* Get big enough correction for S error to let the particle Part be fitted to XYZ point
312   
313   Double_t d[3] = { XYZ[0]-Part[0], XYZ[1]-Part[1], XYZ[2]-Part[2] };
314   Double_t p2 = Part[3]*Part[3]+Part[4]*Part[4]+Part[5]*Part[5];
315   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.;
316   return sigmaS;
317 }
318
319 void AliKFParticleBase::GetMeasurement( const Double_t XYZ[], Double_t m[], Double_t V[] ) const
320 {
321   //* Get additional covariances V used during measurement
322
323   Double_t b[3];
324   GetFieldValue( XYZ, b );
325   const Double_t kCLight =  0.000299792458;
326   b[0]*=kCLight; b[1]*=kCLight; b[2]*=kCLight;
327
328   Transport( GetDStoPoint(XYZ), m, V );
329
330   Double_t sigmaS = GetSCorrection( m, XYZ );
331
332   Double_t h[6];
333
334   h[0] = m[3]*sigmaS;
335   h[1] = m[4]*sigmaS;
336   h[2] = m[5]*sigmaS;
337   h[3] = ( h[1]*b[2]-h[2]*b[1] )*GetQ();
338   h[4] = ( h[2]*b[0]-h[0]*b[2] )*GetQ();
339   h[5] = ( h[0]*b[1]-h[1]*b[0] )*GetQ();
340     
341   V[ 0]+= h[0]*h[0];
342   V[ 1]+= h[1]*h[0];
343   V[ 2]+= h[1]*h[1];
344   V[ 3]+= h[2]*h[0];
345   V[ 4]+= h[2]*h[1];
346   V[ 5]+= h[2]*h[2];
347
348   V[ 6]+= h[3]*h[0];
349   V[ 7]+= h[3]*h[1];
350   V[ 8]+= h[3]*h[2];
351   V[ 9]+= h[3]*h[3];
352
353   V[10]+= h[4]*h[0];
354   V[11]+= h[4]*h[1];
355   V[12]+= h[4]*h[2];
356   V[13]+= h[4]*h[3];
357   V[14]+= h[4]*h[4];
358
359   V[15]+= h[5]*h[0];
360   V[16]+= h[5]*h[1];
361   V[17]+= h[5]*h[2];
362   V[18]+= h[5]*h[3];
363   V[19]+= h[5]*h[4];
364   V[20]+= h[5]*h[5];
365 }
366
367   
368 void AliKFParticleBase::AddDaughter( const AliKFParticleBase &Daughter )
369 {
370   //* Add daughter 
371
372   if( fNDF<-1 ){ // first daughter -> just copy
373     fNDF   = -1;
374     fQ     =  Daughter.GetQ();
375     for( Int_t i=0; i<7; i++) fP[i] = Daughter.fP[i];
376     for( Int_t i=0; i<28; i++) fC[i] = Daughter.fC[i];
377     fSFromDecay = 0;
378     return;
379   }
380
381   TransportToDecayVertex();
382
383   Double_t b[3]; 
384   Int_t maxIter = 1;
385
386   if( !fIsLinearized ){
387     if( fNDF==-1 ){
388       Double_t ds, ds1;
389       GetDStoParticle(Daughter, ds, ds1);      
390       TransportToDS( ds );
391       Double_t m[8];
392       Double_t mCd[36];       
393       Daughter.Transport( ds1, m, mCd );    
394       fVtxGuess[0] = .5*( fP[0] + m[0] );
395       fVtxGuess[1] = .5*( fP[1] + m[1] );
396       fVtxGuess[2] = .5*( fP[2] + m[2] );
397     } else {
398       fVtxGuess[0] = fP[0];
399       fVtxGuess[1] = fP[1];
400       fVtxGuess[2] = fP[2]; 
401     }
402     maxIter = 3;
403   }
404
405   for( Int_t iter=0; iter<maxIter; iter++ ){
406
407     {
408       GetFieldValue( fVtxGuess, b );
409       const Double_t kCLight =  0.000299792458;
410       b[0]*=kCLight; b[1]*=kCLight; b[2]*=kCLight;
411     }
412
413     Double_t *ffP = fP, *ffC = fC, tmpP[8], tmpC[36];
414     if( fNDF==-1 ){            
415       GetMeasurement( fVtxGuess, tmpP, tmpC );
416       ffP = tmpP;
417       ffC = tmpC;
418     }
419
420     Double_t m[8], mV[36];
421
422     if( Daughter.fC[35]>0 ){
423       Daughter.GetMeasurement( fVtxGuess, m, mV );
424     } else {
425       for( Int_t i=0; i<8; i++ ) m[i] = Daughter.fP[i];
426       for( Int_t i=0; i<36; i++ ) mV[i] = Daughter.fC[i];
427     }
428  
429     //*
430
431     Double_t mS[6];
432     {
433       Double_t mSi[6] = { ffC[0]+mV[0], 
434                           ffC[1]+mV[1], ffC[2]+mV[2], 
435                           ffC[3]+mV[3], ffC[4]+mV[4], ffC[5]+mV[5] };
436      
437       mS[0] = mSi[2]*mSi[5] - mSi[4]*mSi[4];
438       mS[1] = mSi[3]*mSi[4] - mSi[1]*mSi[5];
439       mS[2] = mSi[0]*mSi[5] - mSi[3]*mSi[3];
440       mS[3] = mSi[1]*mSi[4] - mSi[2]*mSi[3];
441       mS[4] = mSi[1]*mSi[3] - mSi[0]*mSi[4];
442       mS[5] = mSi[0]*mSi[2] - mSi[1]*mSi[1];     
443       
444       Double_t s = ( mSi[0]*mS[0] + mSi[1]*mS[1] + mSi[3]*mS[3] );      
445
446       s = ( s > 1.E-20 )  ?1./s :0;       
447       mS[0]*=s;
448       mS[1]*=s;
449       mS[2]*=s;
450       mS[3]*=s;
451       mS[4]*=s;
452       mS[5]*=s;
453     }
454     
455     //* Residual (measured - estimated)
456     
457     Double_t zeta[3] = { m[0]-ffP[0], m[1]-ffP[1], m[2]-ffP[2] };    
458
459     
460     //* CHt = CH' - D'
461     
462     Double_t mCHt0[7], mCHt1[7], mCHt2[7];
463     
464     mCHt0[0]=ffC[ 0] ;       mCHt1[0]=ffC[ 1] ;       mCHt2[0]=ffC[ 3] ;
465     mCHt0[1]=ffC[ 1] ;       mCHt1[1]=ffC[ 2] ;       mCHt2[1]=ffC[ 4] ;
466     mCHt0[2]=ffC[ 3] ;       mCHt1[2]=ffC[ 4] ;       mCHt2[2]=ffC[ 5] ;
467     mCHt0[3]=ffC[ 6]-mV[ 6]; mCHt1[3]=ffC[ 7]-mV[ 7]; mCHt2[3]=ffC[ 8]-mV[ 8];
468     mCHt0[4]=ffC[10]-mV[10]; mCHt1[4]=ffC[11]-mV[11]; mCHt2[4]=ffC[12]-mV[12];
469     mCHt0[5]=ffC[15]-mV[15]; mCHt1[5]=ffC[16]-mV[16]; mCHt2[5]=ffC[17]-mV[17];
470     mCHt0[6]=ffC[21]-mV[21]; mCHt1[6]=ffC[22]-mV[22]; mCHt2[6]=ffC[23]-mV[23];
471   
472     //* Kalman gain K = mCH'*S
473     
474     Double_t k0[7], k1[7], k2[7];
475     
476     for(Int_t i=0;i<7;++i){
477       k0[i] = mCHt0[i]*mS[0] + mCHt1[i]*mS[1] + mCHt2[i]*mS[3];
478       k1[i] = mCHt0[i]*mS[1] + mCHt1[i]*mS[2] + mCHt2[i]*mS[4];
479       k2[i] = mCHt0[i]*mS[3] + mCHt1[i]*mS[4] + mCHt2[i]*mS[5];
480     }
481
482    //* New estimation of the vertex position 
483
484     if( iter<maxIter-1 ){
485       for(Int_t i=0; i<3; ++i) 
486         fVtxGuess[i]= ffP[i] + k0[i]*zeta[0]+k1[i]*zeta[1]+k2[i]*zeta[2];
487       continue;
488     }
489
490     // last itearation -> update the particle
491
492     //* Add the daughter momentum to the particle momentum
493     
494     ffP[ 3] += m[ 3];
495     ffP[ 4] += m[ 4];
496     ffP[ 5] += m[ 5];
497     ffP[ 6] += m[ 6];
498   
499     ffC[ 9] += mV[ 9];
500     ffC[13] += mV[13];
501     ffC[14] += mV[14];
502     ffC[18] += mV[18];
503     ffC[19] += mV[19];
504     ffC[20] += mV[20];
505     ffC[24] += mV[24];
506     ffC[25] += mV[25];
507     ffC[26] += mV[26];
508     ffC[27] += mV[27];
509     
510  
511    //* New estimation of the vertex position r += K*zeta
512     
513     for(Int_t i=0;i<7;++i) 
514       fP[i] = ffP[i] + k0[i]*zeta[0] + k1[i]*zeta[1] + k2[i]*zeta[2];
515     
516     //* New covariance matrix C -= K*(mCH')'
517
518     for(Int_t i=0, k=0;i<7;++i){
519       for(Int_t j=0;j<=i;++j,++k){
520         fC[k] = ffC[k] - (k0[i]*mCHt0[j] + k1[i]*mCHt1[j] + k2[i]*mCHt2[j] );
521       }
522     }
523   
524     //* Calculate Chi^2 
525
526     fNDF  += 2;
527     fQ    +=  Daughter.GetQ();
528     fSFromDecay = 0;    
529     fChi2 += (mS[0]*zeta[0] + mS[1]*zeta[1] + mS[3]*zeta[2])*zeta[0]
530       +      (mS[1]*zeta[0] + mS[2]*zeta[1] + mS[4]*zeta[2])*zeta[1]
531       +      (mS[3]*zeta[0] + mS[4]*zeta[1] + mS[5]*zeta[2])*zeta[2];     
532
533   }
534 }
535
536
537 void AliKFParticleBase::SetProductionVertex( const AliKFParticleBase &Vtx )
538 {
539   //* Set production vertex for the particle, when the particle was not used in the vertex fit
540
541   const Double_t *m = Vtx.fP, *mV = Vtx.fC;
542
543   Bool_t noS = ( fC[35]<=0 ); // no decay length allowed
544
545   if( noS ){ 
546     TransportToDecayVertex();
547     fP[7] = 0;
548     fC[28] = fC[29] = fC[30] = fC[31] = fC[32] = fC[33] = fC[35] = fC[35] = 0;
549   } else {
550     TransportToDS( GetDStoPoint( m ) );
551     fP[7] = -fSFromDecay;
552     Convert(1);
553   }
554
555   Double_t mAi[6];
556
557   InvertSym3( fC, mAi );
558
559   Double_t mB[5][3];
560
561   mB[0][0] = fC[ 6]*mAi[0] + fC[ 7]*mAi[1] + fC[ 8]*mAi[3];
562   mB[0][1] = fC[ 6]*mAi[1] + fC[ 7]*mAi[2] + fC[ 8]*mAi[4];
563   mB[0][2] = fC[ 6]*mAi[3] + fC[ 7]*mAi[4] + fC[ 8]*mAi[5];
564
565   mB[1][0] = fC[10]*mAi[0] + fC[11]*mAi[1] + fC[12]*mAi[3];
566   mB[1][1] = fC[10]*mAi[1] + fC[11]*mAi[2] + fC[12]*mAi[4];
567   mB[1][2] = fC[10]*mAi[3] + fC[11]*mAi[4] + fC[12]*mAi[5];
568
569   mB[2][0] = fC[15]*mAi[0] + fC[16]*mAi[1] + fC[17]*mAi[3];
570   mB[2][1] = fC[15]*mAi[1] + fC[16]*mAi[2] + fC[17]*mAi[4];
571   mB[2][2] = fC[15]*mAi[3] + fC[16]*mAi[4] + fC[17]*mAi[5];
572
573   mB[3][0] = fC[21]*mAi[0] + fC[22]*mAi[1] + fC[23]*mAi[3];
574   mB[3][1] = fC[21]*mAi[1] + fC[22]*mAi[2] + fC[23]*mAi[4];
575   mB[3][2] = fC[21]*mAi[3] + fC[22]*mAi[4] + fC[23]*mAi[5];
576
577   mB[4][0] = fC[28]*mAi[0] + fC[29]*mAi[1] + fC[30]*mAi[3];
578   mB[4][1] = fC[28]*mAi[1] + fC[29]*mAi[2] + fC[30]*mAi[4];
579   mB[4][2] = fC[28]*mAi[3] + fC[29]*mAi[4] + fC[30]*mAi[5];
580
581   Double_t z[3] = { m[0]-fP[0], m[1]-fP[1], m[2]-fP[2] };
582
583   {
584     Double_t mAVi[6] = { fC[0]-mV[0], fC[1]-mV[1], fC[2]-mV[2], 
585                         fC[3]-mV[3], fC[4]-mV[4], fC[5]-mV[5] };
586     
587     if( !InvertSym3( mAVi, mAVi ) ){
588
589       Double_t dChi2 = ( +(mAVi[0]*z[0] + mAVi[1]*z[1] + mAVi[3]*z[2])*z[0]
590                          +(mAVi[1]*z[0] + mAVi[2]*z[1] + mAVi[4]*z[2])*z[1]
591                          +(mAVi[3]*z[0] + mAVi[4]*z[1] + mAVi[5]*z[2])*z[2] );
592       
593       // Take Abs(dChi2) here. Negative value of 'det' or 'dChi2' shows that the particle 
594       // was not used in the production vertex fit
595       
596       fChi2+= TMath::Abs( dChi2 );
597     }
598     fNDF  += 2;
599   }
600   
601   fP[0] = m[0];
602   fP[1] = m[1];
603   fP[2] = m[2];
604   fP[3]+= mB[0][0]*z[0] + mB[0][1]*z[1] + mB[0][2]*z[2];
605   fP[4]+= mB[1][0]*z[0] + mB[1][1]*z[1] + mB[1][2]*z[2];
606   fP[5]+= mB[2][0]*z[0] + mB[2][1]*z[1] + mB[2][2]*z[2];
607   fP[6]+= mB[3][0]*z[0] + mB[3][1]*z[1] + mB[3][2]*z[2];
608   fP[7]+= mB[4][0]*z[0] + mB[4][1]*z[1] + mB[4][2]*z[2];
609   
610   Double_t d0, d1, d2;
611
612   fC[0] = mV[0];
613   fC[1] = mV[1];
614   fC[2] = mV[2];
615   fC[3] = mV[3];
616   fC[4] = mV[4];
617   fC[5] = mV[5];
618
619   d0= mB[0][0]*mV[0] + mB[0][1]*mV[1] + mB[0][2]*mV[3] - fC[ 6];
620   d1= mB[0][0]*mV[1] + mB[0][1]*mV[2] + mB[0][2]*mV[4] - fC[ 7];
621   d2= mB[0][0]*mV[3] + mB[0][1]*mV[4] + mB[0][2]*mV[5] - fC[ 8];
622
623   fC[ 6]+= d0;
624   fC[ 7]+= d1;
625   fC[ 8]+= d2;
626   fC[ 9]+= d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
627
628   d0= mB[1][0]*mV[0] + mB[1][1]*mV[1] + mB[1][2]*mV[3] - fC[10];
629   d1= mB[1][0]*mV[1] + mB[1][1]*mV[2] + mB[1][2]*mV[4] - fC[11];
630   d2= mB[1][0]*mV[3] + mB[1][1]*mV[4] + mB[1][2]*mV[5] - fC[12];
631
632   fC[10]+= d0;
633   fC[11]+= d1;
634   fC[12]+= d2;
635   fC[13]+= d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
636   fC[14]+= d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
637
638   d0= mB[2][0]*mV[0] + mB[2][1]*mV[1] + mB[2][2]*mV[3] - fC[15];
639   d1= mB[2][0]*mV[1] + mB[2][1]*mV[2] + mB[2][2]*mV[4] - fC[16];
640   d2= mB[2][0]*mV[3] + mB[2][1]*mV[4] + mB[2][2]*mV[5] - fC[17];
641
642   fC[15]+= d0;
643   fC[16]+= d1;
644   fC[17]+= d2;
645   fC[18]+= d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
646   fC[19]+= d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
647   fC[20]+= d0*mB[2][0] + d1*mB[2][1] + d2*mB[2][2];
648
649   d0= mB[3][0]*mV[0] + mB[3][1]*mV[1] + mB[3][2]*mV[3] - fC[21];
650   d1= mB[3][0]*mV[1] + mB[3][1]*mV[2] + mB[3][2]*mV[4] - fC[22];
651   d2= mB[3][0]*mV[3] + mB[3][1]*mV[4] + mB[3][2]*mV[5] - fC[23];
652
653   fC[21]+= d0;
654   fC[22]+= d1;
655   fC[23]+= d2;
656   fC[24]+= d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
657   fC[25]+= d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
658   fC[26]+= d0*mB[2][0] + d1*mB[2][1] + d2*mB[2][2];
659   fC[27]+= d0*mB[3][0] + d1*mB[3][1] + d2*mB[3][2];
660
661   d0= mB[4][0]*mV[0] + mB[4][1]*mV[1] + mB[4][2]*mV[3] - fC[28];
662   d1= mB[4][0]*mV[1] + mB[4][1]*mV[2] + mB[4][2]*mV[4] - fC[29];
663   d2= mB[4][0]*mV[3] + mB[4][1]*mV[4] + mB[4][2]*mV[5] - fC[30];
664
665   fC[28]+= d0;
666   fC[29]+= d1;
667   fC[30]+= d2;
668   fC[31]+= d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
669   fC[32]+= d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
670   fC[33]+= d0*mB[2][0] + d1*mB[2][1] + d2*mB[2][2];
671   fC[34]+= d0*mB[3][0] + d1*mB[3][1] + d2*mB[3][2];
672   fC[35]+= d0*mB[4][0] + d1*mB[4][1] + d2*mB[4][2];
673   
674   if( noS ){ 
675     fP[7] = 0;
676     fC[28] = fC[29] = fC[30] = fC[31] = fC[32] = fC[33] = fC[35] = fC[35] = 0;
677   } else {
678     TransportToDS( fP[7] );
679     Convert(0);
680   }
681
682   fSFromDecay = 0;
683 }
684
685
686
687 void AliKFParticleBase::SetMassConstraint( Double_t Mass, Double_t SigmaMass )
688 {  
689   //* Set hard( SigmaMass=0 ) or soft (SigmaMass>0) mass constraint 
690
691   Double_t m2 = Mass*Mass;            // measurement, weighted by Mass 
692   Double_t s2 = m2*SigmaMass*SigmaMass; // sigma^2
693
694   Double_t p2 = fP[3]*fP[3] + fP[4]*fP[4] + fP[5]*fP[5]; 
695   Double_t e0 = TMath::Sqrt(m2+p2);
696
697   Double_t mH[8];
698   mH[0] = mH[1] = mH[2] = 0.;
699   mH[3] = -2*fP[3]; 
700   mH[4] = -2*fP[4]; 
701   mH[5] = -2*fP[5]; 
702   mH[6] =  2*fP[6];//e0;
703   mH[7] = 0; 
704
705   Double_t zeta = e0*e0 - e0*fP[6];
706   zeta = m2 - (fP[6]*fP[6]-p2);
707   
708   Double_t mCHt[8], s2_est=0;
709   for( Int_t i=0; i<8; ++i ){
710     mCHt[i] = 0.0;
711     for (Int_t j=0;j<8;++j) mCHt[i] += Cij(i,j)*mH[j];
712     s2_est += mH[i]*mCHt[i];
713   }
714   
715   if( s2_est<1.e-20 ) return; // calculated mass error is already 0, 
716                               // the particle can not be constrained on mass
717
718   Double_t w2 = 1./( s2 + s2_est );
719   fChi2 += zeta*zeta*w2;
720   fNDF  += 1;
721   for( Int_t i=0, ii=0; i<8; ++i ){
722     Double_t ki = mCHt[i]*w2;
723     fP[i]+= ki*zeta;
724     for(Int_t j=0;j<=i;++j) fC[ii++] -= ki*mCHt[j];    
725   }
726 }
727
728
729 void AliKFParticleBase::SetNoDecayLength()
730 {  
731   //* Set no decay length for resonances
732
733   TransportToDecayVertex();
734
735   Double_t h[8];
736   h[0] = h[1] = h[2] = h[3] = h[4] = h[5] = h[6] = 0;
737   h[7] = 1; 
738
739   Double_t zeta = 0 - fP[7];
740   for(Int_t i=0;i<8;++i) zeta -= h[i]*(fP[i]-fP[i]);
741   
742   Double_t s = fC[35];   
743   if( s>1.e-20 ){
744     s = 1./s;
745     fChi2 += zeta*zeta*s;
746     fNDF  += 1;
747     for( Int_t i=0, ii=0; i<7; ++i ){
748       Double_t ki = fC[28+i]*s;
749       fP[i]+= ki*zeta;
750       for(Int_t j=0;j<=i;++j) fC[ii++] -= ki*fC[28+j];    
751     }
752   }
753   fP[7] = 0;
754   fC[28] = fC[29] = fC[30] = fC[31] = fC[32] = fC[33] = fC[35] = fC[35] = 0;
755 }
756
757
758 void AliKFParticleBase::Construct( const AliKFParticleBase* vDaughters[], Int_t NDaughters,
759                                    const AliKFParticleBase *Parent,  Double_t Mass, Bool_t IsConstrained         )
760
761   //* Full reconstruction in one go
762
763   Int_t maxIter = 1;
764   bool wasLinearized = fIsLinearized;
765   if( !fIsLinearized || IsConstrained ){
766     //fVtxGuess[0] = fVtxGuess[1] = fVtxGuess[2] = 0;  //!!!!
767     fVtxGuess[0] = GetX();
768     fVtxGuess[1] = GetY();
769     fVtxGuess[2] = GetZ();
770     fIsLinearized = 1;
771     maxIter = 3;
772   }
773
774   Double_t constraintC[6];
775
776   if( IsConstrained ){
777     for(Int_t i=0;i<6;++i) constraintC[i]=fC[i];
778   } else {
779     for(Int_t i=0;i<6;++i) constraintC[i]=0.;
780     constraintC[0] = constraintC[2] = constraintC[5] = 100.;    
781   }
782
783
784   for( Int_t iter=0; iter<maxIter; iter++ ){
785     fAtProductionVertex = 0;
786     fSFromDecay = 0;
787     fP[0] = fVtxGuess[0];
788     fP[1] = fVtxGuess[1];
789     fP[2] = fVtxGuess[2];
790     fP[3] = 0;
791     fP[4] = 0;
792     fP[5] = 0;
793     fP[6] = 0;
794     fP[7] = 0;
795    
796     for(Int_t i=0;i<6; ++i) fC[i]=constraintC[i];
797     for(Int_t i=6;i<36;++i) fC[i]=0.;
798     fC[35] = 1.;
799     
800     fNDF  = IsConstrained ?0 :-3;
801     fChi2 =  0.;
802     fQ = 0;
803
804     for( Int_t itr =0; itr<NDaughters; itr++ ){
805       AddDaughter( *vDaughters[itr] );    
806     }
807     if( iter<maxIter-1){
808       for( Int_t i=0; i<3; i++ ) fVtxGuess[i] = fP[i];  
809     }
810   }
811   fIsLinearized = wasLinearized;    
812
813   if( Mass>=0 ) SetMassConstraint( Mass );
814   if( Parent ) SetProductionVertex( *Parent );
815 }
816
817
818 void AliKFParticleBase::Convert( bool ToProduction )
819 {
820   //* Tricky function - convert the particle error along its trajectory to 
821   //* the value which corresponds to its production/decay vertex
822   //* It is done by combination of the error of decay length with the position errors
823
824   Double_t fld[3]; 
825   {
826     GetFieldValue( fP, fld );
827     const Double_t kCLight =  fQ*0.000299792458;
828     fld[0]*=kCLight; fld[1]*=kCLight; fld[2]*=kCLight;
829   }
830
831   Double_t h[6];
832   
833   h[0] = fP[3];
834   h[1] = fP[4];
835   h[2] = fP[5];
836   if( ToProduction ){ h[0]=-h[0]; h[1]=-h[1]; h[2]=-h[2]; } 
837   h[3] = h[1]*fld[2]-h[2]*fld[1];
838   h[4] = h[2]*fld[0]-h[0]*fld[2];
839   h[5] = h[0]*fld[1]-h[1]*fld[0];
840   
841   Double_t c;
842
843   c = fC[28]+h[0]*fC[35];
844   fC[ 0]+= h[0]*(c+fC[28]);
845   fC[28] = c;
846
847   fC[ 1]+= h[1]*fC[28] + h[0]*fC[29];
848   c = fC[29]+h[1]*fC[35];
849   fC[ 2]+= h[1]*(c+fC[29]);
850   fC[29] = c;
851
852   fC[ 3]+= h[2]*fC[28] + h[0]*fC[30];
853   fC[ 4]+= h[2]*fC[29] + h[1]*fC[30];
854   c = fC[30]+h[2]*fC[35];
855   fC[ 5]+= h[2]*(c+fC[30]);
856   fC[30] = c;
857
858   fC[ 6]+= h[3]*fC[28] + h[0]*fC[31];
859   fC[ 7]+= h[3]*fC[29] + h[1]*fC[31];
860   fC[ 8]+= h[3]*fC[30] + h[2]*fC[31];
861   c = fC[31]+h[3]*fC[35];
862   fC[ 9]+= h[3]*(c+fC[31]);
863   fC[31] = c;
864   
865   fC[10]+= h[4]*fC[28] + h[0]*fC[32];
866   fC[11]+= h[4]*fC[29] + h[1]*fC[32];
867   fC[12]+= h[4]*fC[30] + h[2]*fC[32];
868   fC[13]+= h[4]*fC[31] + h[3]*fC[32];
869   c = fC[32]+h[4]*fC[35];
870   fC[14]+= h[4]*(c+fC[32]);
871   fC[32] = c;
872   
873   fC[15]+= h[5]*fC[28] + h[0]*fC[33];
874   fC[16]+= h[5]*fC[29] + h[1]*fC[33];
875   fC[17]+= h[5]*fC[30] + h[2]*fC[33];
876   fC[18]+= h[5]*fC[31] + h[3]*fC[33];
877   fC[19]+= h[5]*fC[32] + h[4]*fC[33];
878   c = fC[33]+h[5]*fC[35];
879   fC[20]+= h[5]*(c+fC[33]);
880   fC[33] = c;
881
882   fC[21]+= h[0]*fC[34];
883   fC[22]+= h[1]*fC[34];
884   fC[23]+= h[2]*fC[34];
885   fC[24]+= h[3]*fC[34];
886   fC[25]+= h[4]*fC[34];
887   fC[26]+= h[5]*fC[34];
888 }
889
890
891 void AliKFParticleBase::TransportToDecayVertex()
892 {
893   //* Transport the particle to its decay vertex 
894
895   if( fSFromDecay != 0 ) TransportToDS( -fSFromDecay );
896   if( fAtProductionVertex ) Convert(0);
897   fAtProductionVertex = 0;
898 }
899
900 void AliKFParticleBase::TransportToProductionVertex()
901 {
902   //* Transport the particle to its production vertex 
903   
904   if( fSFromDecay != -fP[7] ) TransportToDS( -fSFromDecay-fP[7] );
905   if( !fAtProductionVertex ) Convert( 1 );
906   fAtProductionVertex = 1;
907 }
908
909
910 void AliKFParticleBase::TransportToDS( Double_t dS )
911
912   //* Transport the particle on dS parameter (SignedPath/Momentum) 
913  
914   Transport( dS, fP, fC );
915   fSFromDecay+= dS;
916 }
917
918
919 Double_t AliKFParticleBase::GetDStoPointLine( const Double_t xyz[] ) const 
920 {
921   //* Get dS to a certain space point without field
922
923   Double_t p2 = fP[3]*fP[3] + fP[4]*fP[4] + fP[5]*fP[5];  
924   if( p2<1.e-4 ) p2 = 1;
925   return ( fP[3]*(xyz[0]-fP[0]) + fP[4]*(xyz[1]-fP[1]) + fP[5]*(xyz[2]-fP[2]) )/p2;
926 }
927
928
929 Double_t AliKFParticleBase::GetDStoPointBz( Double_t B, const Double_t xyz[] ) 
930   const
931
932   
933   //* Get dS to a certain space point for Bz field
934   const Double_t kCLight = 0.000299792458;
935   Double_t bq = B*fQ*kCLight;
936   Double_t pt2 = fP[3]*fP[3] + fP[4]*fP[4];
937   if( pt2<1.e-4 ) return 0;
938   Double_t dx = xyz[0] - fP[0];
939   Double_t dy = xyz[1] - fP[1]; 
940   Double_t a = dx*fP[3]+dy*fP[4];
941   Double_t dS;
942
943   if( TMath::Abs(bq)<1.e-8 ) dS = a/pt2;  
944   else dS =  TMath::ATan2( bq*a, pt2 + bq*(dy*fP[3] -dx*fP[4]) )/bq;
945
946   if(0){
947
948     Double_t px = fP[3];
949     Double_t py = fP[4];
950     Double_t pz = fP[5];
951     Double_t ss[2], g[2][5];
952   
953     ss[0] = dS;
954     ss[1] = -dS;
955     for( Int_t i=0; i<2; i++){
956       Double_t bs = bq*ss[i];
957       Double_t c = TMath::Cos(bs), s = TMath::Sin(bs);
958       Double_t cB,sB;
959       if( TMath::Abs(bq)>1.e-8){
960         cB= (1-c)/bq;     
961         sB= s/bq;  
962       }else{
963         const Double_t kOvSqr6 = 1./TMath::Sqrt(6.);
964         sB = (1.-bs*kOvSqr6)*(1.+bs*kOvSqr6)*ss[i];
965         cB = .5*sB*bs;
966       }
967       g[i][0] = fP[0] + sB*px + cB*py;
968       g[i][1] = fP[1] - cB*px + sB*py;
969       g[i][2] = fP[2] + ss[i]*pz;
970       g[i][3] =       + c*px + s*py;
971       g[i][4] =       - s*px + c*py;      
972     }
973
974     Int_t i=0;
975   
976     Double_t dMin = 1.e10;
977     for( Int_t j=0; j<2; j++){
978       Double_t xx = g[j][0]-xyz[0];
979       Double_t yy = g[j][1]-xyz[1];
980       Double_t zz = g[j][2]-xyz[2];
981       Double_t d = xx*xx + yy*yy + zz*zz;
982       if( d<dMin ){
983         dMin = d;
984         i = j;
985       }
986     }
987
988     dS = ss[i];
989
990     Double_t x= g[i][0], y= g[i][1], z= g[i][2], ppx= g[i][3], ppy= g[i][4];      
991     Double_t ddx = x-xyz[0];
992     Double_t ddy = y-xyz[1];
993     Double_t ddz = z-xyz[2];
994     Double_t c = ddx*ppx  + ddy*ppy  + ddz*pz ;
995     Double_t pp2 = ppx*ppx + ppy*ppy + pz*pz;    
996     if( TMath::Abs(pp2)>1.e-8 ){
997       dS+=c/pp2;
998     }
999   }
1000   return dS;
1001 }
1002
1003
1004 void AliKFParticleBase::GetDStoParticleBz( Double_t B, const AliKFParticleBase &p, 
1005                                            Double_t &DS, Double_t &DS1 ) 
1006   const
1007
1008   //* Get dS to another particle for Bz field
1009   Double_t px = fP[3];
1010   Double_t py = fP[4];
1011   Double_t pz = fP[5];
1012
1013   Double_t px1 = p.fP[3];
1014   Double_t py1 = p.fP[4];
1015   Double_t pz1 = p.fP[5];
1016
1017   const Double_t kCLight = 0.000299792458;
1018
1019   Double_t bq = B*fQ*kCLight;
1020   Double_t bq1 = B*p.fQ*kCLight;
1021   Double_t s=0, ds=0, s1=0, ds1=0;
1022   
1023   if( TMath::Abs(bq)>1.e-8 || TMath::Abs(bq1)>1.e-8 ){
1024
1025     Double_t dx = (p.fP[0] - fP[0]);
1026     Double_t dy = (p.fP[1] - fP[1]);
1027     Double_t d2 = (dx*dx+dy*dy);
1028     
1029     Double_t p2  = (px *px  + py *py); 
1030     Double_t p21 = (px1*px1 + py1*py1);
1031     
1032     Double_t a = (px*py1 - py*px1);
1033     Double_t b = (px*px1 + py*py1);
1034     
1035     Double_t ldx = bq*bq1*dx - bq1*py + bq*py1 ;
1036     Double_t ldy = bq*bq1*dy + bq1*px - bq*px1 ;
1037     Double_t l2 = ldx*ldx + ldy*ldy;
1038     
1039     Double_t cS = bq1*p2 + bq*bq1*(dy* px - dx* py) -  bq*b;
1040     Double_t cS1= bq*p21 - bq*bq1*(dy*px1 - dx*py1) - bq1*b;
1041
1042     Double_t ca  = bq*bq*bq1*d2  +2*( cS + bq*bq*(py1*dx-px1*dy)) ;
1043     Double_t ca1 = bq*bq1*bq1*d2 +2*( cS1 - bq1*bq1*(py*dx-px*dy)) ;  
1044   
1045     Double_t sa = 4*l2*p2 - ca*ca;
1046     Double_t sa1 = 4*l2*p21 - ca1*ca1;
1047
1048     if(sa<0) sa=0;
1049     if(sa1<0)sa1=0;
1050
1051     if( TMath::Abs(bq)>1.e-8){
1052       s  = TMath::ATan2(   bq*( bq1*(dx*px +dy*py) + a ) , cS )/bq;
1053       ds = TMath::ATan2(TMath::Sqrt(sa),ca)/bq;
1054     } else {
1055       s = ( (dx*px + dy*py) + (py*px1-px*py1)/bq1)/p2;
1056       ds = s*s - (d2-2*(px1*dy-py1*dx)/bq1)/p2; 
1057       if( ds<0 ) ds = 0;
1058       ds = TMath::Sqrt(ds);   
1059     }
1060     
1061     if( TMath::Abs(bq1)>1.e-8){
1062       s1 = TMath::ATan2( -bq1*( bq*(dx*px1+dy*py1) + a), cS1 )/bq1;
1063       ds1 = TMath::ATan2(TMath::Sqrt(sa1),ca1)/bq1;  
1064     } else {
1065       s1 = (-(dx*px1 + dy*py1) + (py*px1-px*py1)/bq)/p21;
1066       ds1 = s1*s1 - (d2+2*(px*dy-py*dx)/bq)/p21; 
1067       if( ds1<0 ) ds1 = 0;
1068       ds1 = TMath::Sqrt(ds1);
1069     }
1070   }
1071
1072   Double_t ss[2], ss1[2], g[2][5],g1[2][5];
1073   
1074   ss[0] = s + ds;
1075   ss[1] = s - ds;
1076   ss1[0] = s1 + ds1;
1077   ss1[1] = s1 - ds1;
1078   for( Int_t i=0; i<2; i++){
1079     Double_t bs = bq*ss[i];
1080     Double_t c = TMath::Cos(bs), sss = TMath::Sin(bs);
1081     Double_t cB,sB;
1082     if( TMath::Abs(bq)>1.e-8){
1083       cB= (1-c)/bq;     
1084       sB= sss/bq;  
1085     }else{
1086       const Double_t kOvSqr6 = 1./TMath::Sqrt(6.);
1087       sB = (1.-bs*kOvSqr6)*(1.+bs*kOvSqr6)*ss[i];
1088       cB = .5*sB*bs;
1089     }
1090     g[i][0] = fP[0] + sB*px + cB*py;
1091     g[i][1] = fP[1] - cB*px + sB*py;
1092     g[i][2] = fP[2] + ss[i]*pz;
1093     g[i][3] =       + c*px + sss*py;
1094     g[i][4] =       - sss*px + c*py;
1095
1096     bs = bq1*ss1[i];  
1097     c =  TMath::Cos(bs); sss = TMath::Sin(bs);
1098     if( TMath::Abs(bq1)>1.e-8){
1099       cB= (1-c)/bq1;   
1100       sB= sss/bq1;  
1101     }else{
1102       const Double_t kOvSqr6 = 1./TMath::Sqrt(6.);
1103       sB = (1.-bs*kOvSqr6)*(1.+bs*kOvSqr6)*ss1[i];
1104       cB = .5*sB*bs;
1105     }
1106       
1107     g1[i][0] = p.fP[0] + sB*px1 + cB*py1;
1108     g1[i][1] = p.fP[1] - cB*px1 + sB*py1;
1109     g1[i][2] = p.fP[2] + ss[i]*pz1;
1110     g1[i][3] =         + c*px1 + sss*py1;
1111     g1[i][4] =         - sss*px1 + c*py1;
1112   }
1113
1114   Int_t i=0, i1=0;
1115   
1116   Double_t dMin = 1.e10;
1117   for( Int_t j=0; j<2; j++){
1118     for( Int_t j1=0; j1<2; j1++){
1119       Double_t xx = g[j][0]-g1[j1][0];
1120       Double_t yy = g[j][1]-g1[j1][1];
1121       Double_t zz = g[j][2]-g1[j1][2];
1122       Double_t d = xx*xx + yy*yy + zz*zz;
1123       if( d<dMin ){
1124         dMin = d;
1125         i = j;
1126         i1 = j1;
1127       }
1128     }
1129   }  
1130
1131   DS = ss[i];
1132   DS1 = ss1[i1];
1133   if(0){
1134     Double_t x= g[i][0], y= g[i][1], z= g[i][2], ppx= g[i][3], ppy= g[i][4];  
1135     Double_t x1=g1[i1][0], y1= g1[i1][1], z1= g1[i1][2], ppx1= g1[i1][3], ppy1= g1[i1][4];  
1136     Double_t dx = x1-x;
1137     Double_t dy = y1-y;
1138     Double_t dz = z1-z;
1139     Double_t a = ppx*ppx1 + ppy*ppy1 + pz*pz1;
1140     Double_t b = dx*ppx1 + dy*ppy1 + dz*pz1;
1141     Double_t c = dx*ppx  + dy*ppy  + dz*pz ;
1142     Double_t pp2 = ppx*ppx + ppy*ppy + pz*pz;
1143     Double_t pp21= ppx1*ppx1 + ppy1*ppy1 + pz1*pz1;
1144     Double_t det = pp2*pp21 - a*a;    
1145     if( TMath::Abs(det)>1.e-8 ){
1146       DS+=(a*b-pp21*c)/det;
1147       DS1+=(a*c-pp2*b)/det;
1148     }
1149   }
1150 }
1151
1152
1153
1154 void AliKFParticleBase::TransportCBM( Double_t dS, 
1155                                  Double_t P[], Double_t C[] ) const
1156 {  
1157   //* Transport the particle on dS, output to P[],C[], for CBM field
1158  
1159   if( fQ==0 ){
1160     TransportLine( dS, P, C );
1161     return;
1162   }
1163
1164   const Double_t kCLight = 0.000299792458;
1165
1166   Double_t c = fQ*kCLight;
1167
1168   // construct coefficients 
1169
1170   Double_t 
1171     px   = fP[3],
1172     py   = fP[4],
1173     pz   = fP[5];
1174       
1175   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;
1176
1177   { // get field integrals
1178
1179     Double_t fld[3][3];   
1180     Double_t p0[3], p1[3], p2[3];
1181
1182     // line track approximation
1183
1184     p0[0] = fP[0];
1185     p0[1] = fP[1];
1186     p0[2] = fP[2];
1187   
1188     p2[0] = fP[0] + px*dS;
1189     p2[1] = fP[1] + py*dS;
1190     p2[2] = fP[2] + pz*dS;
1191   
1192     p1[0] = 0.5*(p0[0]+p2[0]);
1193     p1[1] = 0.5*(p0[1]+p2[1]);
1194     p1[2] = 0.5*(p0[2]+p2[2]);
1195
1196     // first order track approximation
1197     {
1198       GetFieldValue( p0, fld[0] );
1199       GetFieldValue( p1, fld[1] );
1200       GetFieldValue( p2, fld[2] );
1201
1202       Double_t ssy1 = ( 7*fld[0][1] + 6*fld[1][1]-fld[2][1] )*c*dS*dS/96.;
1203       Double_t ssy2 = (   fld[0][1] + 2*fld[1][1]         )*c*dS*dS/6.;
1204
1205       p1[0] -= ssy1*pz;
1206       p1[2] += ssy1*px;
1207       p2[0] -= ssy2*pz;
1208       p2[2] += ssy2*px;   
1209     }
1210
1211     GetFieldValue( p0, fld[0] );
1212     GetFieldValue( p1, fld[1] );
1213     GetFieldValue( p2, fld[2] );
1214     
1215     sx = c*( fld[0][0] + 4*fld[1][0] + fld[2][0] )*dS/6.;
1216     sy = c*( fld[0][1] + 4*fld[1][1] + fld[2][1] )*dS/6.;
1217     sz = c*( fld[0][2] + 4*fld[1][2] + fld[2][2] )*dS/6.;
1218
1219     ssx = c*( fld[0][0] + 2*fld[1][0])*dS*dS/6.;
1220     ssy = c*( fld[0][1] + 2*fld[1][1])*dS*dS/6.;
1221     ssz = c*( fld[0][2] + 2*fld[1][2])*dS*dS/6.;
1222
1223     Double_t c2[3][3]    =   { {  5, -4, -1},{  44,  80,  -4},{ 11, 44, 5} }; // /=360.    
1224     Double_t cc2[3][3]    =   { { 38,  8, -4},{ 148, 208, -20},{  3, 36, 3} }; // /=2520.
1225     for(Int_t n=0; n<3; n++)
1226       for(Int_t m=0; m<3; m++) 
1227         {
1228           syz += c2[n][m]*fld[n][1]*fld[m][2];
1229           ssyz += cc2[n][m]*fld[n][1]*fld[m][2];
1230         }
1231  
1232     syz  *= c*c*dS*dS/360.;
1233     ssyz  *= c*c*dS*dS*dS/2520.;
1234     
1235     syy  = c*( fld[0][1] + 4*fld[1][1] + fld[2][1] )*dS;
1236     syyy = syy*syy*syy / 1296;
1237     syy  = syy*syy/72;
1238
1239     ssyy = ( fld[0][1]*( 38*fld[0][1] + 156*fld[1][1]  -   fld[2][1] )+
1240             fld[1][1]*(              208*fld[1][1]  +16*fld[2][1] )+
1241             fld[2][1]*(                             3*fld[2][1] )  
1242             )*dS*dS*dS*c*c/2520.;
1243     ssyyy = 
1244       (
1245        fld[0][1]*( fld[0][1]*( 85*fld[0][1] + 526*fld[1][1]  - 7*fld[2][1] )+
1246                  fld[1][1]*(             1376*fld[1][1]  +84*fld[2][1] )+
1247                  fld[2][1]*(                            19*fld[2][1] )  )+
1248        fld[1][1]*( fld[1][1]*(             1376*fld[1][1] +256*fld[2][1] )+
1249                  fld[2][1]*(                            62*fld[2][1] )  )+
1250        fld[2][1]*fld[2][1]  *(                             3*fld[2][1] )       
1251        )*dS*dS*dS*dS*c*c*c/90720.;    
1252  
1253   }
1254
1255   Double_t mJ[8][8];
1256   for( Int_t i=0; i<8; i++ ) for( Int_t j=0; j<8; j++) mJ[i][j]=0;
1257
1258   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;
1259   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;
1260   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;
1261   
1262   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;
1263   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;
1264   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;
1265   mJ[6][6] = mJ[7][7] = 1;
1266   
1267   P[0] = fP[0] + mJ[0][3]*px + mJ[0][4]*py + mJ[0][5]*pz;
1268   P[1] = fP[1] + mJ[1][3]*px + mJ[1][4]*py + mJ[1][5]*pz;
1269   P[2] = fP[2] + mJ[2][3]*px + mJ[2][4]*py + mJ[2][5]*pz;
1270   P[3] =        mJ[3][3]*px + mJ[3][4]*py + mJ[3][5]*pz;
1271   P[4] =        mJ[4][3]*px + mJ[4][4]*py + mJ[4][5]*pz;
1272   P[5] =        mJ[5][3]*px + mJ[5][4]*py + mJ[5][5]*pz;
1273   P[6] = fP[6];
1274   P[7] = fP[7];
1275
1276   MultQSQt( mJ[0], fC, C);
1277
1278 }
1279
1280
1281 void AliKFParticleBase::TransportBz( Double_t b, Double_t t,
1282                                      Double_t p[], Double_t e[] ) const 
1283
1284   //* Transport the particle on dS, output to P[],C[], for Bz field
1285  
1286   const Double_t kCLight = 0.000299792458;
1287   b = b*fQ*kCLight;
1288   Double_t bs= b*t;
1289   Double_t s = TMath::Sin(bs), c = TMath::Cos(bs);
1290   Double_t sB, cB;
1291   if( TMath::Abs(bs)>1.e-10){
1292     sB= s/b;
1293     cB= (1-c)/b;
1294   }else{
1295     const Double_t kOvSqr6 = 1./TMath::Sqrt(6.);
1296     sB = (1.-bs*kOvSqr6)*(1.+bs*kOvSqr6)*t;
1297     cB = .5*sB*bs;
1298   }
1299   
1300   Double_t px = fP[3];
1301   Double_t py = fP[4];
1302   Double_t pz = fP[5];
1303   
1304   p[0] = fP[0] + sB*px + cB*py;
1305   p[1] = fP[1] - cB*px + sB*py;
1306   p[2] = fP[2] +  t*pz;
1307   p[3] =          c*px + s*py;
1308   p[4] =         -s*px + c*py;
1309   p[5] = fP[5];
1310   p[6] = fP[6];
1311   p[7] = fP[7];
1312
1313   /* 
1314   Double_t mJ[8][8] = { {1,0,0,   sB, cB,  0, 0, 0 },
1315                         {0,1,0,  -cB, sB,  0, 0, 0 },
1316                         {0,0,1,    0,  0,  t, 0, 0 },
1317                         {0,0,0,    c,  s,  0, 0, 0 },
1318                         {0,0,0,   -s,  c,  0, 0, 0 },
1319                         {0,0,0,    0,  0,  1, 0, 0 },
1320                         {0,0,0,    0,  0,  0, 1, 0 },
1321                         {0,0,0,    0,  0,  0, 0, 1 }  };
1322   Double_t mA[8][8];
1323   for( Int_t k=0,i=0; i<8; i++)
1324     for( Int_t j=0; j<=i; j++, k++ ) mA[i][j] = mA[j][i] = fC[k]; 
1325
1326   Double_t mJC[8][8];
1327   for( Int_t i=0; i<8; i++ )
1328     for( Int_t j=0; j<8; j++ ){
1329       mJC[i][j]=0;
1330       for( Int_t k=0; k<8; k++ ) mJC[i][j]+=mJ[i][k]*mA[k][j];
1331     }
1332   
1333   for( Int_t k=0,i=0; i<8; i++)
1334     for( Int_t j=0; j<=i; j++, k++ ){
1335       e[k] = 0;
1336       for( Int_t l=0; l<8; l++ ) e[k]+=mJC[i][l]*mJ[j][l];
1337     }
1338   
1339   return;
1340   */
1341
1342   Double_t 
1343     c6=fC[6], c7=fC[7], c8=fC[8], c17=fC[17], c18=fC[18],
1344     c24 = fC[24], c31 = fC[31];
1345
1346   Double_t 
1347     cBC13 = cB*fC[13],
1348     mJC13 = c7 - cB*fC[9] + sB*fC[13],
1349     mJC14 = fC[11] - cBC13 + sB*fC[14],
1350     mJC23 = c8 + t*c18,
1351     mJC24 = fC[12] + t*fC[19],
1352     mJC33 = c*fC[9] + s*fC[13],
1353     mJC34 = c*fC[13] + s*fC[14],
1354     mJC43 = -s*fC[9] + c*fC[13],
1355     mJC44 = -s*fC[13] + c*fC[14];
1356
1357
1358   e[0]= fC[0] + 2*(sB*c6 + cB*fC[10]) + (sB*fC[9] + 2*cBC13)*sB + cB*cB*fC[14];
1359   e[1]= fC[1] - cB*c6 + sB*fC[10] + mJC13*sB + mJC14*cB;
1360   e[2]= fC[2] - cB*c7 + sB*fC[11] - mJC13*cB + mJC14*sB;
1361   e[3]= fC[3] + t*fC[15] + mJC23*sB + mJC24*cB;
1362   e[4]= fC[4] + t*fC[16] - mJC23*cB + mJC24*sB;
1363
1364   e[15]= fC[15] + c18*sB + fC[19]*cB;
1365   e[16]= fC[16] - c18*cB + fC[19]*sB;
1366   e[17]= c17 + fC[20]*t;
1367   e[18]= c18*c + fC[19]*s;
1368   e[19]= -c18*s + fC[19]*c;
1369
1370   e[5]= fC[5] + (c17 + e[17] )*t;
1371
1372   e[6]= c*c6 + s*fC[10] + mJC33*sB + mJC34*cB;
1373   e[7]= c*c7 + s*fC[11] - mJC33*cB + mJC34*sB;
1374   e[8]= c*c8 + s*fC[12] + e[18]*t;
1375   e[9]= mJC33*c + mJC34*s;
1376   e[10]= -s*c6 + c*fC[10] + mJC43*sB + mJC44*cB;
1377
1378     
1379   e[11]= -s*c7 + c*fC[11] - mJC43*cB + mJC44*sB;
1380   e[12]= -s*c8 + c*fC[12] + e[19]*t;
1381   e[13]= mJC43*c + mJC44*s;
1382   e[14]= -mJC43*s + mJC44*c;
1383   e[20]= fC[20];
1384   e[21]= fC[21] + fC[25]*cB + c24*sB;
1385   e[22]= fC[22] - c24*cB + fC[25]*sB;
1386   e[23]= fC[23] + fC[26]*t;
1387   e[24]= c*c24 + s*fC[25];
1388   e[25]= c*fC[25] - c24*s;
1389   e[26]= fC[26];
1390   e[27]= fC[27];
1391   e[28]= fC[28] + fC[32]*cB + c31*sB;
1392   e[29]= fC[29] - c31*cB + fC[32]*sB;
1393   e[30]= fC[30] + fC[33]*t;
1394   e[31]= c*c31 + s*fC[32];
1395   e[32]= c*fC[32] - s*c31;
1396   e[33]= fC[33];
1397   e[34]= fC[34];
1398   e[35]= fC[35];     
1399 }
1400
1401
1402 Double_t AliKFParticleBase::GetDistanceFromVertex( const AliKFParticleBase &Vtx ) const
1403 {
1404   //* Calculate distance from vertex [cm]
1405
1406   return GetDistanceFromVertex( Vtx.fP );
1407 }
1408
1409 Double_t AliKFParticleBase::GetDistanceFromVertex( const Double_t vtx[] ) const
1410 {
1411   //* Calculate distance from vertex [cm]
1412
1413   Double_t mP[8], mC[36];  
1414   Transport( GetDStoPoint(vtx), mP, mC );
1415   Double_t d[3]={ vtx[0]-mP[0], vtx[1]-mP[1], vtx[2]-mP[2]};
1416   return TMath::Sqrt( d[0]*d[0]+d[1]*d[1]+d[2]*d[2] );
1417 }
1418
1419 Double_t AliKFParticleBase::GetDistanceFromParticle( const AliKFParticleBase &p ) 
1420   const
1421
1422   //* Calculate distance to other particle [cm]
1423
1424   Double_t dS, dS1;
1425   GetDStoParticle( p, dS, dS1 );   
1426   Double_t mP[8], mC[36], mP1[8], mC1[36];
1427   Transport( dS, mP, mC ); 
1428   p.Transport( dS1, mP1, mC1 ); 
1429   Double_t dx = mP[0]-mP1[0]; 
1430   Double_t dy = mP[1]-mP1[1]; 
1431   Double_t dz = mP[2]-mP1[2]; 
1432   return TMath::Sqrt(dx*dx+dy*dy+dz*dz);
1433 }
1434
1435 Double_t AliKFParticleBase::GetDeviationFromVertex( const AliKFParticleBase &Vtx ) const
1436 {
1437   //* Calculate sqrt(Chi2/ndf) deviation from vertex
1438
1439   return GetDeviationFromVertex( Vtx.fP, Vtx.fC );
1440 }
1441
1442
1443 Double_t AliKFParticleBase::GetDeviationFromVertex( const Double_t v[], const Double_t Cv[] ) const
1444 {
1445   //* Calculate sqrt(Chi2/ndf) deviation from vertex
1446   //* v = [xyz], Cv=[Cxx,Cxy,Cyy,Cxz,Cyz,Czz]-covariance matrix
1447
1448   Double_t mP[8];
1449   Double_t mC[36];
1450   
1451   Transport( GetDStoPoint(v), mP, mC );  
1452
1453   Double_t d[3]={ v[0]-mP[0], v[1]-mP[1], v[2]-mP[2]};
1454
1455   Double_t sigmaS = .1+10.*TMath::Sqrt( (d[0]*d[0]+d[1]*d[1]+d[2]*d[2])/
1456                                  (mP[3]*mP[3]+mP[4]*mP[4]+mP[5]*mP[5])  );
1457
1458    
1459   Double_t h[3] = { mP[3]*sigmaS, mP[4]*sigmaS, mP[5]*sigmaS };       
1460   
1461   Double_t mSi[6] = 
1462     { mC[0] +h[0]*h[0], 
1463       mC[1] +h[1]*h[0], mC[2] +h[1]*h[1], 
1464       mC[3] +h[2]*h[0], mC[4] +h[2]*h[1], mC[5] +h[2]*h[2] };
1465
1466   if( Cv ){
1467     mSi[0]+=Cv[0];
1468     mSi[1]+=Cv[1];
1469     mSi[2]+=Cv[2];
1470     mSi[3]+=Cv[3];
1471     mSi[4]+=Cv[4];
1472     mSi[5]+=Cv[5];
1473   }
1474   
1475   Double_t mS[6];
1476
1477   mS[0] = mSi[2]*mSi[5] - mSi[4]*mSi[4];
1478   mS[1] = mSi[3]*mSi[4] - mSi[1]*mSi[5];
1479   mS[2] = mSi[0]*mSi[5] - mSi[3]*mSi[3];
1480   mS[3] = mSi[1]*mSi[4] - mSi[2]*mSi[3];
1481   mS[4] = mSi[1]*mSi[3] - mSi[0]*mSi[4];
1482   mS[5] = mSi[0]*mSi[2] - mSi[1]*mSi[1];         
1483       
1484   Double_t s = ( mSi[0]*mS[0] + mSi[1]*mS[1] + mSi[3]*mS[3] );
1485   s = ( s > 1.E-20 )  ?1./s :0;   
1486
1487   return TMath::Sqrt( TMath::Abs(s*( ( mS[0]*d[0] + mS[1]*d[1] + mS[3]*d[2])*d[0]
1488                    +(mS[1]*d[0] + mS[2]*d[1] + mS[4]*d[2])*d[1]
1489                    +(mS[3]*d[0] + mS[4]*d[1] + mS[5]*d[2])*d[2] ))/2);
1490 }
1491
1492
1493 Double_t AliKFParticleBase::GetDeviationFromParticle( const AliKFParticleBase &p ) 
1494   const
1495
1496   //* Calculate sqrt(Chi2/ndf) deviation from other particle
1497
1498   Double_t dS, dS1;
1499   GetDStoParticle( p, dS, dS1 );   
1500   Double_t mP1[8], mC1[36];
1501   p.Transport( dS1, mP1, mC1 ); 
1502
1503   Double_t d[3]={ fP[0]-mP1[0], fP[1]-mP1[1], fP[2]-mP1[2]};
1504
1505   Double_t sigmaS = .1+10.*TMath::Sqrt( (d[0]*d[0]+d[1]*d[1]+d[2]*d[2])/
1506                                         (mP1[3]*mP1[3]+mP1[4]*mP1[4]+mP1[5]*mP1[5])  );
1507
1508   Double_t h[3] = { mP1[3]*sigmaS, mP1[4]*sigmaS, mP1[5]*sigmaS };       
1509   
1510   mC1[0] +=h[0]*h[0];
1511   mC1[1] +=h[1]*h[0]; 
1512   mC1[2] +=h[1]*h[1]; 
1513   mC1[3] +=h[2]*h[0]; 
1514   mC1[4] +=h[2]*h[1];
1515   mC1[5] +=h[2]*h[2];
1516
1517   return GetDeviationFromVertex( mP1, mC1 )*TMath::Sqrt(2./1.);
1518 }
1519
1520
1521
1522 void AliKFParticleBase::SubtractFromVertex(  AliKFParticleBase &Vtx ) const
1523 {    
1524   //* Subtract the particle from the vertex  
1525
1526   Double_t fld[3];  
1527   {
1528     GetFieldValue( Vtx.fP, fld );
1529     const Double_t kCLight =  0.000299792458;
1530     fld[0]*=kCLight; fld[1]*=kCLight; fld[2]*=kCLight;
1531   }
1532
1533   Double_t m[8];
1534   Double_t mCm[36];
1535
1536   if( Vtx.fIsLinearized ){
1537     GetMeasurement( Vtx.fVtxGuess, m, mCm );
1538   } else {
1539     GetMeasurement( Vtx.fP, m, mCm );
1540   }
1541   
1542   Double_t mV[6];
1543     
1544   mV[ 0] = mCm[ 0];
1545   mV[ 1] = mCm[ 1];
1546   mV[ 2] = mCm[ 2];
1547   mV[ 3] = mCm[ 3];
1548   mV[ 4] = mCm[ 4];
1549   mV[ 5] = mCm[ 5];
1550      
1551   //* 
1552             
1553   Double_t mS[6];
1554   {
1555     Double_t mSi[6] = { mV[0]-Vtx.fC[0], 
1556                         mV[1]-Vtx.fC[1], mV[2]-Vtx.fC[2], 
1557                         mV[3]-Vtx.fC[3], mV[4]-Vtx.fC[4], mV[5]-Vtx.fC[5] };
1558     
1559     mS[0] = mSi[2]*mSi[5] - mSi[4]*mSi[4];
1560     mS[1] = mSi[3]*mSi[4] - mSi[1]*mSi[5];
1561     mS[2] = mSi[0]*mSi[5] - mSi[3]*mSi[3];
1562     mS[3] = mSi[1]*mSi[4] - mSi[2]*mSi[3];
1563     mS[4] = mSi[1]*mSi[3] - mSi[0]*mSi[4];
1564     mS[5] = mSi[0]*mSi[2] - mSi[1]*mSi[1];       
1565     
1566     Double_t s = ( mSi[0]*mS[0] + mSi[1]*mS[1] + mSi[3]*mS[3] );
1567     s = ( s > 1.E-20 )  ?1./s :0;         
1568     mS[0]*=s;
1569     mS[1]*=s;
1570     mS[2]*=s;
1571     mS[3]*=s;
1572     mS[4]*=s;
1573     mS[5]*=s;
1574   }
1575     
1576   //* Residual (measured - estimated)
1577     
1578   Double_t zeta[3] = { m[0]-Vtx.fP[0], m[1]-Vtx.fP[1], m[2]-Vtx.fP[2] };
1579         
1580   //* mCHt = mCH' - D'
1581     
1582   Double_t mCHt0[3], mCHt1[3], mCHt2[3];
1583     
1584   mCHt0[0]=Vtx.fC[ 0] ;      mCHt1[0]=Vtx.fC[ 1] ;      mCHt2[0]=Vtx.fC[ 3] ;
1585   mCHt0[1]=Vtx.fC[ 1] ;      mCHt1[1]=Vtx.fC[ 2] ;      mCHt2[1]=Vtx.fC[ 4] ;
1586   mCHt0[2]=Vtx.fC[ 3] ;      mCHt1[2]=Vtx.fC[ 4] ;      mCHt2[2]=Vtx.fC[ 5] ;
1587   
1588   //* Kalman gain K = mCH'*S
1589     
1590   Double_t k0[3], k1[3], k2[3];
1591     
1592   for(Int_t i=0;i<3;++i){
1593     k0[i] = mCHt0[i]*mS[0] + mCHt1[i]*mS[1] + mCHt2[i]*mS[3];
1594     k1[i] = mCHt0[i]*mS[1] + mCHt1[i]*mS[2] + mCHt2[i]*mS[4];
1595     k2[i] = mCHt0[i]*mS[3] + mCHt1[i]*mS[4] + mCHt2[i]*mS[5];
1596   }
1597     
1598   //* New estimation of the vertex position r += K*zeta
1599     
1600   Double_t dChi2 = -(mS[0]*zeta[0] + mS[1]*zeta[1] + mS[3]*zeta[2])*zeta[0]
1601     +      (mS[1]*zeta[0] + mS[2]*zeta[1] + mS[4]*zeta[2])*zeta[1]
1602     +      (mS[3]*zeta[0] + mS[4]*zeta[1] + mS[5]*zeta[2])*zeta[2];
1603
1604   if( Vtx.fChi2 - dChi2 < 0 ) return;
1605
1606   for(Int_t i=0;i<3;++i) 
1607     Vtx.fP[i] -= k0[i]*zeta[0] + k1[i]*zeta[1] + k2[i]*zeta[2];       
1608     
1609   //* New covariance matrix C -= K*(mCH')'
1610     
1611   for(Int_t i=0, k=0;i<3;++i){
1612     for(Int_t j=0;j<=i;++j,++k) 
1613       Vtx.fC[k] += k0[i]*mCHt0[j] + k1[i]*mCHt1[j] + k2[i]*mCHt2[j];
1614   }
1615     
1616   //* Calculate Chi^2 
1617
1618   Vtx.fNDF  -= 2;
1619   Vtx.fChi2 -= dChi2;
1620 }
1621
1622
1623
1624 void AliKFParticleBase::TransportLine( Double_t dS, 
1625                                        Double_t P[], Double_t C[] ) const 
1626 {
1627   //* Transport the particle as a straight line
1628
1629   P[0] = fP[0] + dS*fP[3];
1630   P[1] = fP[1] + dS*fP[4];
1631   P[2] = fP[2] + dS*fP[5];
1632   P[3] = fP[3];
1633   P[4] = fP[4];
1634   P[5] = fP[5];
1635   P[6] = fP[6];
1636   P[7] = fP[7];
1637  
1638   Double_t c6  = fC[ 6] + dS*fC[ 9];
1639   Double_t c11 = fC[11] + dS*fC[14];
1640   Double_t c17 = fC[17] + dS*fC[20];
1641   Double_t sc13 = dS*fC[13];
1642   Double_t sc18 = dS*fC[18];
1643   Double_t sc19 = dS*fC[19];
1644
1645   C[ 0] = fC[ 0] + dS*( fC[ 6] + c6  );
1646   C[ 2] = fC[ 2] + dS*( fC[11] + c11 );
1647   C[ 5] = fC[ 5] + dS*( fC[17] + c17 );
1648
1649   C[ 7] = fC[ 7] + sc13;
1650   C[ 8] = fC[ 8] + sc18;
1651   C[ 9] = fC[ 9];
1652
1653   C[12] = fC[12] + sc19;
1654
1655   C[ 1] = fC[ 1] + dS*( fC[10] + C[ 7] );
1656   C[ 3] = fC[ 3] + dS*( fC[15] + C[ 8] );
1657   C[ 4] = fC[ 4] + dS*( fC[16] + C[12] ); 
1658   C[ 6] = c6;
1659
1660   C[10] = fC[10] + sc13;
1661   C[11] = c11;
1662
1663   C[13] = fC[13];
1664   C[14] = fC[14];
1665   C[15] = fC[15] + sc18;
1666   C[16] = fC[16] + sc19;
1667   C[17] = c17;
1668   
1669   C[18] = fC[18];
1670   C[19] = fC[19];
1671   C[20] = fC[20];
1672   C[21] = fC[21] + dS*fC[24];
1673   C[22] = fC[22] + dS*fC[25];
1674   C[23] = fC[23] + dS*fC[26];
1675
1676   C[24] = fC[24];
1677   C[25] = fC[25];
1678   C[26] = fC[26];
1679   C[27] = fC[27];
1680   C[28] = fC[28] + dS*fC[31];
1681   C[29] = fC[29] + dS*fC[32];
1682   C[30] = fC[30] + dS*fC[33];
1683
1684   C[31] = fC[31];
1685   C[32] = fC[32];
1686   C[33] = fC[33];
1687   C[34] = fC[34];
1688   C[35] = fC[35]; 
1689 }
1690
1691
1692 void AliKFParticleBase::ConstructGammaBz( const AliKFParticleBase &daughter1,
1693                                           const AliKFParticleBase &daughter2, double Bz  )
1694
1695   //* Create gamma
1696   
1697   const AliKFParticleBase *daughters[2] = { &daughter1, &daughter2};
1698
1699   double v0[3];
1700   
1701   if( !fIsLinearized ){
1702     Double_t ds, ds1;
1703     Double_t m[8];
1704     Double_t mCd[36];       
1705     daughter1.GetDStoParticle(daughter2, ds, ds1);      
1706     daughter1.Transport( ds, m, mCd );
1707     v0[0] = m[0];
1708     v0[1] = m[1];
1709     v0[2] = m[2];
1710     daughter2.Transport( ds1, m, mCd );
1711     v0[0] = .5*( v0[0] + m[0] );
1712     v0[1] = .5*( v0[1] + m[1] );
1713     v0[2] = .5*( v0[2] + m[2] );
1714   } else {
1715     v0[0] = fVtxGuess[0];
1716     v0[1] = fVtxGuess[1];
1717     v0[2] = fVtxGuess[2];
1718   }
1719
1720   fAtProductionVertex = 0;
1721   fSFromDecay = 0;
1722   fP[0] = v0[0];
1723   fP[1] = v0[1];
1724   fP[2] = v0[2];
1725   fP[3] = 0;
1726   fP[4] = 0;
1727   fP[5] = 0;
1728   fP[6] = 0;
1729   fP[7] = 0;
1730
1731   
1732   // fit daughters to the vertex guess
1733   
1734   double daughterP[2][8], daughterC[2][36];
1735   double vtxMom[2][3];
1736
1737   {  
1738     for( int id=0; id<2; id++ ){
1739     
1740       double *p = daughterP[id];
1741       double *mC = daughterC[id];
1742       
1743       daughters[id]->GetMeasurement( v0, p, mC );
1744       
1745       Double_t mAi[6];
1746       InvertSym3(mC, mAi );
1747       
1748       Double_t mB[3][3];
1749       
1750       mB[0][0] = mC[ 6]*mAi[0] + mC[ 7]*mAi[1] + mC[ 8]*mAi[3];
1751       mB[0][1] = mC[ 6]*mAi[1] + mC[ 7]*mAi[2] + mC[ 8]*mAi[4];
1752       mB[0][2] = mC[ 6]*mAi[3] + mC[ 7]*mAi[4] + mC[ 8]*mAi[5];
1753       
1754       mB[1][0] = mC[10]*mAi[0] + mC[11]*mAi[1] + mC[12]*mAi[3];
1755       mB[1][1] = mC[10]*mAi[1] + mC[11]*mAi[2] + mC[12]*mAi[4];
1756       mB[1][2] = mC[10]*mAi[3] + mC[11]*mAi[4] + mC[12]*mAi[5];
1757       
1758       mB[2][0] = mC[15]*mAi[0] + mC[16]*mAi[1] + mC[17]*mAi[3];
1759       mB[2][1] = mC[15]*mAi[1] + mC[16]*mAi[2] + mC[17]*mAi[4];
1760       mB[2][2] = mC[15]*mAi[3] + mC[16]*mAi[4] + mC[17]*mAi[5];
1761     
1762       Double_t z[3] = { v0[0]-p[0], v0[1]-p[1], v0[2]-p[2] };
1763     
1764       vtxMom[id][0] = p[3] + mB[0][0]*z[0] + mB[0][1]*z[1] + mB[0][2]*z[2];
1765       vtxMom[id][1] = p[4] + mB[1][0]*z[0] + mB[1][1]*z[1] + mB[1][2]*z[2];
1766       vtxMom[id][2] = p[5] + mB[2][0]*z[0] + mB[2][1]*z[1] + mB[2][2]*z[2];
1767                   
1768       daughters[id]->Transport( daughters[id]->GetDStoPoint(v0), p, mC );
1769       
1770     }      
1771
1772   } // fit daughters to guess
1773
1774
1775     // fit new vertex
1776   {
1777       
1778     double mpx0 =  vtxMom[0][0]+vtxMom[1][0];
1779     double mpy0 =  vtxMom[0][1]+vtxMom[1][1];
1780     double mpt0 = TMath::Sqrt(mpx0*mpx0 + mpy0*mpy0);
1781     // double a0 = TMath::ATan2(mpy0,mpx0);
1782     
1783     double ca0 = mpx0/mpt0;
1784     double sa0 = mpy0/mpt0;
1785     double r[3] = { v0[0], v0[1], v0[2] };
1786     double mC[3][3] = {{10000., 0 ,   0  },
1787                        {0,  10000.,   0  },
1788                        {0,     0, 10000. } };
1789     double chi2=0;
1790     
1791     for( int id=0; id<2; id++ ){                
1792       const Double_t kCLight = 0.000299792458;
1793       Double_t q = Bz*daughters[id]->GetQ()*kCLight;
1794       Double_t px0 = vtxMom[id][0];
1795       Double_t py0 = vtxMom[id][1];
1796       Double_t pz0 = vtxMom[id][2];
1797       Double_t pt0 = TMath::Sqrt(px0*px0+py0*py0);
1798       Double_t mG[3][6], mB[3], mH[3][3];
1799       // r = {vx,vy,vz};
1800       // m = {x,y,z,Px,Py,Pz};
1801       // V = daughter.C
1802       // G*m + B = H*r;
1803       // q*x + Py - q*vx - sin(a)*Pt = 0
1804       // q*y - Px - q*vy + cos(a)*Pt = 0
1805       // (Px*cos(a) + Py*sin(a) ) (vz -z) - Pz( cos(a)*(vx-x) + sin(a)*(vy-y)) = 0
1806       
1807       mG[0][0] = q;
1808       mG[0][1] = 0;
1809       mG[0][2] = 0;
1810       mG[0][3] =   -sa0*px0/pt0;
1811       mG[0][4] = 1 -sa0*py0/pt0;
1812       mG[0][5] = 0;     
1813       mH[0][0] = q;
1814       mH[0][1] = 0;
1815       mH[0][2] = 0;      
1816       mB[0] = py0 - sa0*pt0 - mG[0][3]*px0 - mG[0][4]*py0 ;
1817       
1818       // q*y - Px - q*vy + cos(a)*Pt = 0
1819       
1820       mG[1][0] = 0;
1821       mG[1][1] = q;
1822       mG[1][2] = 0;
1823       mG[1][3] = -1 + ca0*px0/pt0;
1824       mG[1][4] =    + ca0*py0/pt0;
1825       mG[1][5] = 0;      
1826       mH[1][0] = 0;
1827       mH[1][1] = q;
1828       mH[1][2] = 0;      
1829       mB[1] = -px0 + ca0*pt0 - mG[1][3]*px0 - mG[1][4]*py0 ;
1830       
1831       // (Px*cos(a) + Py*sin(a) ) (z -vz) - Pz( cos(a)*(x-vx) + sin(a)*(y-vy)) = 0
1832       
1833       mG[2][0] = -pz0*ca0;
1834       mG[2][1] = -pz0*sa0;
1835       mG[2][2] =  px0*ca0 + py0*sa0;
1836       mG[2][3] = 0;
1837       mG[2][4] = 0;
1838       mG[2][5] = 0;
1839       
1840       mH[2][0] = mG[2][0];
1841       mH[2][1] = mG[2][1];
1842       mH[2][2] = mG[2][2];
1843       
1844       mB[2] = 0;
1845       
1846       // fit the vertex
1847
1848       // V = GVGt
1849
1850       double mGV[3][6];
1851       double mV[6];
1852       double m[3];
1853       for( int i=0; i<3; i++ ){
1854         m[i] = mB[i];
1855         for( int k=0; k<6; k++ ) m[i]+=mG[i][k]*daughterP[id][k];
1856       }
1857       for( int i=0; i<3; i++ ){
1858         for( int j=0; j<6; j++ ){
1859           mGV[i][j] = 0;
1860           for( int k=0; k<6; k++ ) mGV[i][j]+=mG[i][k]*daughterC[id][ IJ(k,j) ];
1861         }
1862       }
1863       for( int i=0, k=0; i<3; i++ ){
1864         for( int j=0; j<=i; j++,k++ ){
1865           mV[k] = 0;
1866           for( int l=0; l<6; l++ ) mV[k]+=mGV[i][l]*mG[j][l];
1867         }
1868       }
1869       
1870       
1871       //* CHt
1872       
1873       Double_t mCHt[3][3];
1874       Double_t mHCHt[6];
1875       Double_t mHr[3];
1876       for( int i=0; i<3; i++ ){   
1877         mHr[i] = 0;
1878         for( int k=0; k<3; k++ ) mHr[i]+= mH[i][k]*r[k];
1879       }
1880       
1881       for( int i=0; i<3; i++ ){
1882         for( int j=0; j<3; j++){
1883           mCHt[i][j] = 0;
1884           for( int k=0; k<3; k++ ) mCHt[i][j]+= mC[i][k]*mH[j][k];
1885         }
1886       }
1887
1888       for( int i=0, k=0; i<3; i++ ){
1889         for( int j=0; j<=i; j++, k++ ){
1890           mHCHt[k] = 0;
1891           for( int l=0; l<3; l++ ) mHCHt[k]+= mH[i][l]*mCHt[l][j];
1892         }
1893       }
1894       
1895       Double_t mS[6] = { mHCHt[0]+mV[0], 
1896                           mHCHt[1]+mV[1], mHCHt[2]+mV[2], 
1897                           mHCHt[3]+mV[3], mHCHt[4]+mV[4], mHCHt[5]+mV[5]    };  
1898       
1899
1900       InvertSym3(mS,mS);
1901         
1902       //* Residual (measured - estimated)
1903     
1904       Double_t zeta[3] = { m[0]-mHr[0], m[1]-mHr[1], m[2]-mHr[2] };
1905             
1906       //* Kalman gain K = mCH'*S
1907     
1908       Double_t k[3][3];
1909       
1910       for(Int_t i=0;i<3;++i){
1911         k[i][0] = mCHt[i][0]*mS[0] + mCHt[i][1]*mS[1] + mCHt[i][2]*mS[3];
1912         k[i][1] = mCHt[i][0]*mS[1] + mCHt[i][1]*mS[2] + mCHt[i][2]*mS[4];
1913         k[i][2] = mCHt[i][0]*mS[3] + mCHt[i][1]*mS[4] + mCHt[i][2]*mS[5];
1914       }
1915
1916       //* New estimation of the vertex position r += K*zeta
1917     
1918       for(Int_t i=0;i<3;++i) 
1919         r[i] = r[i] + k[i][0]*zeta[0] + k[i][1]*zeta[1] + k[i][2]*zeta[2];
1920       
1921       //* New covariance matrix C -= K*(mCH')'
1922
1923       for(Int_t i=0;i<3;++i){
1924         for(Int_t j=0;j<=i;++j){
1925           mC[i][j] = mC[i][j] - (k[i][0]*mCHt[j][0] + k[i][1]*mCHt[j][1] + k[i][2]*mCHt[j][2]);
1926           mC[j][i] = mC[i][j];
1927         }
1928       }
1929
1930         
1931       //* Calculate Chi^2 
1932         
1933       chi2 += ( ( mS[0]*zeta[0] + mS[1]*zeta[1] + mS[3]*zeta[2] )*zeta[0]
1934                 +(mS[1]*zeta[0] + mS[2]*zeta[1] + mS[4]*zeta[2] )*zeta[1]
1935                 +(mS[3]*zeta[0] + mS[4]*zeta[1] + mS[5]*zeta[2] )*zeta[2]  );  
1936     }
1937     
1938     // store vertex
1939     
1940     fNDF  = 2;
1941     fChi2 = chi2;
1942     for( int i=0; i<3; i++ ) fP[i] = r[i];
1943     for( int i=0,k=0; i<3; i++ ){
1944       for( int j=0; j<=i; j++,k++ ){
1945         fC[k] = mC[i][j];
1946       }
1947     }
1948   }
1949
1950   // now fit daughters to the vertex
1951   
1952   fQ     =  0;
1953   fSFromDecay = 0;    
1954
1955   for(Int_t i=3;i<8;++i) fP[i]=0.;
1956   for(Int_t i=6;i<36;++i) fC[i]=0.;
1957
1958
1959   for( int id=0; id<2; id++ ){
1960
1961     double *p = daughterP[id];
1962     double *mC = daughterC[id];      
1963     daughters[id]->GetMeasurement( v0, p, mC );
1964
1965     const Double_t *m = fP, *mV = fC;
1966     
1967     Double_t mAi[6];
1968     InvertSym3(mC, mAi );
1969
1970     Double_t mB[4][3];
1971
1972     mB[0][0] = mC[ 6]*mAi[0] + mC[ 7]*mAi[1] + mC[ 8]*mAi[3];
1973     mB[0][1] = mC[ 6]*mAi[1] + mC[ 7]*mAi[2] + mC[ 8]*mAi[4];
1974     mB[0][2] = mC[ 6]*mAi[3] + mC[ 7]*mAi[4] + mC[ 8]*mAi[5];
1975     
1976     mB[1][0] = mC[10]*mAi[0] + mC[11]*mAi[1] + mC[12]*mAi[3];
1977     mB[1][1] = mC[10]*mAi[1] + mC[11]*mAi[2] + mC[12]*mAi[4];
1978     mB[1][2] = mC[10]*mAi[3] + mC[11]*mAi[4] + mC[12]*mAi[5];
1979     
1980     mB[2][0] = mC[15]*mAi[0] + mC[16]*mAi[1] + mC[17]*mAi[3];
1981     mB[2][1] = mC[15]*mAi[1] + mC[16]*mAi[2] + mC[17]*mAi[4];
1982     mB[2][2] = mC[15]*mAi[3] + mC[16]*mAi[4] + mC[17]*mAi[5];
1983     
1984     mB[3][0] = mC[21]*mAi[0] + mC[22]*mAi[1] + mC[23]*mAi[3];
1985     mB[3][1] = mC[21]*mAi[1] + mC[22]*mAi[2] + mC[23]*mAi[4];
1986     mB[3][2] = mC[21]*mAi[3] + mC[22]*mAi[4] + mC[23]*mAi[5];    
1987
1988
1989     Double_t z[3] = { m[0]-p[0], m[1]-p[1], m[2]-p[2] };
1990
1991     {
1992       Double_t mAV[6] = { mC[0]-mV[0], mC[1]-mV[1], mC[2]-mV[2], 
1993                           mC[3]-mV[3], mC[4]-mV[4], mC[5]-mV[5] };
1994       
1995       Double_t mAVi[6];
1996       if( !InvertSym3(mAV, mAVi) ){
1997         Double_t dChi2 = ( +(mAVi[0]*z[0] + mAVi[1]*z[1] + mAVi[3]*z[2])*z[0]
1998                            +(mAVi[1]*z[0] + mAVi[2]*z[1] + mAVi[4]*z[2])*z[1]
1999                            +(mAVi[3]*z[0] + mAVi[4]*z[1] + mAVi[5]*z[2])*z[2] );
2000         fChi2+= TMath::Abs( dChi2 );
2001       }
2002       fNDF  += 2;
2003     }
2004
2005     //* Add the daughter momentum to the particle momentum
2006  
2007     fP[3]+= p[3] + mB[0][0]*z[0] + mB[0][1]*z[1] + mB[0][2]*z[2];
2008     fP[4]+= p[4] + mB[1][0]*z[0] + mB[1][1]*z[1] + mB[1][2]*z[2];
2009     fP[5]+= p[5] + mB[2][0]*z[0] + mB[2][1]*z[1] + mB[2][2]*z[2];
2010     fP[6]+= p[6] + mB[3][0]*z[0] + mB[3][1]*z[1] + mB[3][2]*z[2];
2011   
2012     Double_t d0, d1, d2;
2013    
2014     d0= mB[0][0]*mV[0] + mB[0][1]*mV[1] + mB[0][2]*mV[3] - mC[ 6];
2015     d1= mB[0][0]*mV[1] + mB[0][1]*mV[2] + mB[0][2]*mV[4] - mC[ 7];
2016     d2= mB[0][0]*mV[3] + mB[0][1]*mV[4] + mB[0][2]*mV[5] - mC[ 8];
2017
2018     //fC[6]+= mC[ 6] + d0;
2019     //fC[7]+= mC[ 7] + d1;
2020     //fC[8]+= mC[ 8] + d2;
2021     fC[9]+= mC[ 9] + d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
2022
2023     d0= mB[1][0]*mV[0] + mB[1][1]*mV[1] + mB[1][2]*mV[3] - mC[10];
2024     d1= mB[1][0]*mV[1] + mB[1][1]*mV[2] + mB[1][2]*mV[4] - mC[11];
2025     d2= mB[1][0]*mV[3] + mB[1][1]*mV[4] + mB[1][2]*mV[5] - mC[12];
2026
2027     //fC[10]+= mC[10]+ d0;
2028     //fC[11]+= mC[11]+ d1;
2029     //fC[12]+= mC[12]+ d2;
2030     fC[13]+= mC[13]+ d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
2031     fC[14]+= mC[14]+ d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
2032
2033     d0= mB[2][0]*mV[0] + mB[2][1]*mV[1] + mB[2][2]*mV[3] - mC[15];
2034     d1= mB[2][0]*mV[1] + mB[2][1]*mV[2] + mB[2][2]*mV[4] - mC[16];
2035     d2= mB[2][0]*mV[3] + mB[2][1]*mV[4] + mB[2][2]*mV[5] - mC[17];
2036
2037     //fC[15]+= mC[15]+ d0;
2038     //fC[16]+= mC[16]+ d1;
2039     //fC[17]+= mC[17]+ d2;
2040     fC[18]+= mC[18]+ d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
2041     fC[19]+= mC[19]+ d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
2042     fC[20]+= mC[20]+ d0*mB[2][0] + d1*mB[2][1] + d2*mB[2][2];
2043
2044     d0= mB[3][0]*mV[0] + mB[3][1]*mV[1] + mB[3][2]*mV[3] - mC[21];
2045     d1= mB[3][0]*mV[1] + mB[3][1]*mV[2] + mB[3][2]*mV[4] - mC[22];
2046     d2= mB[3][0]*mV[3] + mB[3][1]*mV[4] + mB[3][2]*mV[5] - mC[23];
2047
2048     //fC[21]+= mC[21] + d0;
2049     //fC[22]+= mC[22] + d1;
2050     //fC[23]+= mC[23] + d2;
2051     fC[24]+= mC[24] + d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
2052     fC[25]+= mC[25] + d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
2053     fC[26]+= mC[26] + d0*mB[2][0] + d1*mB[2][1] + d2*mB[2][2];
2054     fC[27]+= mC[27] + d0*mB[3][0] + d1*mB[3][1] + d2*mB[3][2];
2055   }
2056
2057   SetMassConstraint(0,0);
2058   
2059 }
2060
2061 Bool_t AliKFParticleBase::InvertSym3( const Double_t A[], Double_t Ai[] )
2062 {
2063   //* Invert symmetric matric stored in low-triagonal form 
2064
2065   bool ret = 0;
2066   double a0 = A[0], a1 = A[1], a2 = A[2], a3 = A[3];
2067
2068   Ai[0] = a2*A[5] - A[4]*A[4];
2069   Ai[1] = a3*A[4] - a1*A[5];
2070   Ai[3] = a1*A[4] - a2*a3;
2071   Double_t det = (a0*Ai[0] + a1*Ai[1] + a3*Ai[3]);
2072   if( det>1.e-20 ) det = 1./det;    
2073   else{ 
2074     det = 0;
2075     ret = 1;
2076   }
2077   Ai[0] *= det;
2078   Ai[1] *= det;
2079   Ai[3] *= det;
2080   Ai[2] = ( a0*A[5] - a3*a3 )*det;
2081   Ai[4] = ( a1*a3 - a0*A[4] )*det;
2082   Ai[5] = ( a0*a2 - a1*a1 )*det;
2083   return ret;
2084 }
2085
2086 void AliKFParticleBase::MultQSQt( const Double_t Q[], const Double_t S[], Double_t SOut[] )
2087 {
2088   //* Matrix multiplication Q*S*Q^T, Q - square matrix, S - symmetric
2089
2090   const Int_t kN= 8;
2091   Double_t mA[kN*kN];
2092   
2093   for( Int_t i=0, ij=0; i<kN; i++ ){
2094     for( Int_t j=0; j<kN; j++, ++ij ){
2095       mA[ij] = 0 ;
2096       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];
2097     }
2098   }
2099     
2100   for( Int_t i=0; i<kN; i++ ){
2101     for( Int_t j=0; j<=i; j++ ){
2102       Int_t ij = ( j<=i ) ? i*(i+1)/2+j :j*(j+1)/2+i;
2103       SOut[ij] = 0 ;
2104       for( Int_t k=0; k<kN; k++ )  SOut[ij] += Q[ i*kN+k ] * mA[ k*kN+j ];
2105     }
2106   }
2107 }
2108
2109
2110 // 72-charachters line to define the printer border
2111 //3456789012345678901234567890123456789012345678901234567890123456789012
2112