]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/ESD/AliKFParticleBase.cxx
fabs -> TMath::Abs
[u/mrichter/AliRoot.git] / STEER / ESD / AliKFParticleBase.cxx
1 //---------------------------------------------------------------------------------
2 // Implementation of the AliKFParticleBase class
3 // .
4 // @author  S.Gorbunov, I.Kisel, I.Kulakov, M.Zyzak
5 // @version 1.0
6 // @since   13.05.07
7 // 
8 // Class to reconstruct and store the decayed particle parameters.
9 // The method is described in CBM-SOFT note 2007-003, 
10 // ``Reconstruction of decayed particles based on the Kalman filter'', 
11 // http://www.gsi.de/documents/DOC-2007-May-14-1.pdf
12 //
13 // This class 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 #include <iostream>
23 ClassImp(AliKFParticleBase)
24
25
26 AliKFParticleBase::AliKFParticleBase() :fQ(0), fNDF(-3), fChi2(0), fSFromDecay(0), fAtProductionVertex(0), fIsLinearized(0),
27                                         fConstructMethod(2), SumDaughterMass(0), fMassHypo(-1)
28
29   //* Constructor 
30
31   Initialize();
32 }
33
34 void AliKFParticleBase::Initialize( const Double_t Param[], const Double_t Cov[], Int_t Charge, Double_t Mass )
35 {
36   // Constructor from "cartesian" track, particle mass hypothesis should be provided
37   //
38   // Param[6] = { X, Y, Z, Px, Py, Pz } - position and momentum
39   // Cov [21] = lower-triangular part of the covariance matrix:
40   //
41   //                (  0  .  .  .  .  . )
42   //                (  1  2  .  .  .  . )
43   //  Cov. matrix = (  3  4  5  .  .  . ) - numbering of covariance elements in Cov[]
44   //                (  6  7  8  9  .  . )
45   //                ( 10 11 12 13 14  . )
46   //                ( 15 16 17 18 19 20 )
47
48
49   for( Int_t i=0; i<6 ; i++ ) fP[i] = Param[i];
50   for( Int_t i=0; i<21; i++ ) fC[i] = Cov[i];
51
52   Double_t energy = TMath::Sqrt( Mass*Mass + fP[3]*fP[3] + fP[4]*fP[4] + fP[5]*fP[5]);
53   fP[6] = energy;
54   fP[7] = 0;
55   fQ = Charge;
56   fNDF = 0;
57   fChi2 = 0;
58   fAtProductionVertex = 0;
59   fIsLinearized = 0;
60   fSFromDecay = 0;
61
62   Double_t energyInv = 1./energy;
63   Double_t 
64     h0 = fP[3]*energyInv,
65     h1 = fP[4]*energyInv,
66     h2 = fP[5]*energyInv;
67
68   fC[21] = h0*fC[ 6] + h1*fC[10] + h2*fC[15];
69   fC[22] = h0*fC[ 7] + h1*fC[11] + h2*fC[16];
70   fC[23] = h0*fC[ 8] + h1*fC[12] + h2*fC[17];
71   fC[24] = h0*fC[ 9] + h1*fC[13] + h2*fC[18];
72   fC[25] = h0*fC[13] + h1*fC[14] + h2*fC[19];
73   fC[26] = h0*fC[18] + h1*fC[19] + h2*fC[20];
74   fC[27] = ( h0*h0*fC[ 9] + h1*h1*fC[14] + h2*h2*fC[20] 
75              + 2*(h0*h1*fC[13] + h0*h2*fC[18] + h1*h2*fC[19] ) );
76   for( Int_t i=28; i<36; i++ ) fC[i] = 0;
77   fC[35] = 1.;
78
79   SumDaughterMass = Mass;
80   fMassHypo = Mass;
81 }
82
83 void AliKFParticleBase::Initialize()
84 {
85   //* Initialise covariance matrix and set current parameters to 0.0 
86
87   for( Int_t i=0; i<8; i++) fP[i] = 0;
88   for(Int_t i=0;i<36;++i) fC[i]=0.;
89   fC[0] = fC[2] = fC[5] = 100.;
90   fC[35] = 1.;
91   fNDF  = -3;
92   fChi2 =  0.;
93   fQ = 0;
94   fSFromDecay = 0;
95   fAtProductionVertex = 0;
96   fVtxGuess[0]=fVtxGuess[1]=fVtxGuess[2]=0.;
97   fIsLinearized = 0;
98   SumDaughterMass = 0;
99   fMassHypo = -1;
100 }
101
102 void AliKFParticleBase::SetVtxGuess( Double_t x, Double_t y, Double_t z )
103 {
104   //* Set decay vertex parameters for linearisation 
105
106   fVtxGuess[0] = x;
107   fVtxGuess[1] = y;
108   fVtxGuess[2] = z;
109   fIsLinearized = 1;
110 }
111
112 Int_t AliKFParticleBase::GetMomentum( Double_t &p, Double_t &error )  const 
113 {
114   //* Calculate particle momentum
115
116   Double_t x = fP[3];
117   Double_t y = fP[4];
118   Double_t z = fP[5];
119   
120   Double_t x2 = x*x;
121   Double_t y2 = y*y;
122   Double_t z2 = z*z;
123   Double_t p2 = x2+y2+z2;
124   p = TMath::Sqrt(p2);
125   
126   error = (x2*fC[9]+y2*fC[14]+z2*fC[20] + 2*(x*y*fC[13]+x*z*fC[18]+y*z*fC[19]) );
127   if( error>1.e-16 && p>1.e-4 ){
128     error = TMath::Sqrt(error)/p;
129     return 0;
130   }
131   error = 1.e8;
132   return 1;
133 }
134
135 Int_t AliKFParticleBase::GetPt( Double_t &pt, Double_t &error )  const 
136 {
137   //* Calculate particle transverse momentum
138
139   Double_t px = fP[3];
140   Double_t py = fP[4];
141   Double_t px2 = px*px;
142   Double_t py2 = py*py;
143   Double_t pt2 = px2+py2;
144   pt = TMath::Sqrt(pt2);
145   error = (px2*fC[9] + py2*fC[14] + 2*px*py*fC[13] );
146   if( error>0 && pt>1.e-4 ){
147     error = TMath::Sqrt(error)/pt;
148     return 0;
149   }
150   error = 1.e10;
151   return 1;
152 }
153
154 Int_t AliKFParticleBase::GetEta( Double_t &eta, Double_t &error )  const 
155 {
156   //* Calculate particle pseudorapidity
157
158   Double_t px = fP[3];
159   Double_t py = fP[4];
160   Double_t pz = fP[5];
161   Double_t pt2 = px*px + py*py;
162   Double_t p2 = pt2 + pz*pz;
163   Double_t p = TMath::Sqrt(p2);
164   Double_t a = p + pz;
165   Double_t b = p - pz;
166   eta = 1.e10;
167   if( b > 1.e-8 ){
168     Double_t c = a/b;
169     if( c>1.e-8 ) eta = 0.5*TMath::Log(c);
170   }
171   Double_t h3 = -px*pz;
172   Double_t h4 = -py*pz;  
173   Double_t pt4 = pt2*pt2;
174   Double_t p2pt4 = p2*pt4;
175   error = (h3*h3*fC[9] + h4*h4*fC[14] + pt4*fC[20] + 2*( h3*(h4*fC[13] + fC[18]*pt2) + pt2*h4*fC[19] ) );
176
177   if( error>0 && p2pt4>1.e-10 ){
178     error = TMath::Sqrt(error/p2pt4);
179     return 0;
180   }
181
182   error = 1.e10;
183   return 1;
184 }
185
186 Int_t AliKFParticleBase::GetPhi( Double_t &phi, Double_t &error )  const 
187 {
188   //* Calculate particle polar angle
189
190   Double_t px = fP[3];
191   Double_t py = fP[4];
192   Double_t px2 = px*px;
193   Double_t py2 = py*py;
194   Double_t pt2 = px2 + py2;
195   phi = TMath::ATan2(py,px);
196   error = (py2*fC[9] + px2*fC[14] - 2*px*py*fC[13] );
197   if( error>0 && pt2>1.e-4 ){
198     error = TMath::Sqrt(error)/pt2;
199     return 0;
200   }
201   error = 1.e10;
202   return 1;
203 }
204
205 Int_t AliKFParticleBase::GetR( Double_t &r, Double_t &error )  const 
206 {
207   //* Calculate distance to the origin
208
209   Double_t x = fP[0];
210   Double_t y = fP[1];
211   Double_t x2 = x*x;
212   Double_t y2 = y*y;
213   r = TMath::Sqrt(x2 + y2);
214   error = (x2*fC[0] + y2*fC[2] - 2*x*y*fC[1] );
215   if( error>0 && r>1.e-4 ){
216     error = TMath::Sqrt(error)/r;
217     return 0;
218   }
219   error = 1.e10;
220   return 1;
221 }
222
223 Int_t AliKFParticleBase::GetMass( Double_t &m, Double_t &error ) const 
224 {
225   //* Calculate particle mass
226   
227   // s = sigma^2 of m2/2
228
229   Double_t s = (  fP[3]*fP[3]*fC[9] + fP[4]*fP[4]*fC[14] + fP[5]*fP[5]*fC[20] 
230                   + fP[6]*fP[6]*fC[27] 
231                 +2*( + fP[3]*fP[4]*fC[13] + fP[5]*(fP[3]*fC[18] + fP[4]*fC[19]) 
232                      - fP[6]*( fP[3]*fC[24] + fP[4]*fC[25] + fP[5]*fC[26] )   )
233                  ); 
234 //   Double_t m2 = TMath::Abs(fP[6]*fP[6] - fP[3]*fP[3] - fP[4]*fP[4] - fP[5]*fP[5]);
235 //   m  = TMath::Sqrt(m2);
236 //   if( m>1.e-10 ){
237 //     if( s>=0 ){
238 //       error = TMath::Sqrt(s)/m;
239 //       return 0;
240 //     }
241 //   }
242 //   error = 1.e20;
243 //   return 1;
244   Double_t m2 = (fP[6]*fP[6] - fP[3]*fP[3] - fP[4]*fP[4] - fP[5]*fP[5]);
245
246   if(m2<0.)
247   {
248     error = 1.e20;
249     m = -TMath::Sqrt(-m2);
250     return 1;
251   }
252
253   m  = TMath::Sqrt(m2);
254   if( m>1.e-6 ){
255     if( s >= 0 ) {
256       error = TMath::Sqrt(s)/m;
257       return 0;
258     }
259   }
260   else {
261     error = 1.e20;
262     return 0;
263   }
264   error = 1.e20;
265
266   return 1;
267 }
268
269
270 Int_t AliKFParticleBase::GetDecayLength( Double_t &l, Double_t &error ) const 
271 {
272   //* Calculate particle decay length [cm]
273
274   Double_t x = fP[3];
275   Double_t y = fP[4];
276   Double_t z = fP[5];
277   Double_t t = fP[7];
278   Double_t x2 = x*x;
279   Double_t y2 = y*y;
280   Double_t z2 = z*z;
281   Double_t p2 = x2+y2+z2;
282   l = t*TMath::Sqrt(p2);
283   if( p2>1.e-4){
284     error = p2*fC[35] + t*t/p2*(x2*fC[9]+y2*fC[14]+z2*fC[20]
285                                 + 2*(x*y*fC[13]+x*z*fC[18]+y*z*fC[19]) )
286       + 2*t*(x*fC[31]+y*fC[32]+z*fC[33]);
287     error = TMath::Sqrt(TMath::Abs(error));
288     return 0;
289   }
290   error = 1.e20;
291   return 1;
292 }
293
294 Int_t AliKFParticleBase::GetDecayLengthXY( Double_t &l, Double_t &error ) const 
295 {
296   //* Calculate particle decay length in XY projection [cm]
297
298   Double_t x = fP[3];
299   Double_t y = fP[4];
300   Double_t t = fP[7];
301   Double_t x2 = x*x;
302   Double_t y2 = y*y;
303   Double_t pt2 = x2+y2;
304   l = t*TMath::Sqrt(pt2);
305   if( pt2>1.e-4){
306     error = pt2*fC[35] + t*t/pt2*(x2*fC[9]+y2*fC[14] + 2*x*y*fC[13] )
307       + 2*t*(x*fC[31]+y*fC[32]);
308     error = TMath::Sqrt(TMath::Abs(error));
309     return 0;
310   }
311   error = 1.e20;
312   return 1;
313 }
314
315
316 Int_t AliKFParticleBase::GetLifeTime( Double_t &tauC, Double_t &error ) const 
317 {
318   //* Calculate particle decay time [s]
319
320   Double_t m, dm;
321   GetMass( m, dm );
322   Double_t cTM = (-fP[3]*fC[31] - fP[4]*fC[32] - fP[5]*fC[33] + fP[6]*fC[34]);
323   tauC = fP[7]*m;
324   error = m*m*fC[35] + 2*fP[7]*cTM + fP[7]*fP[7]*dm*dm;
325   if( error > 0 ){
326     error = TMath::Sqrt( error );
327     return 0;
328   }
329   error = 1.e20;
330   return 1;
331 }
332
333
334 void AliKFParticleBase::operator +=( const AliKFParticleBase &Daughter )
335 {
336   //* Add daughter via operator+=
337
338   AddDaughter( Daughter );
339 }
340   
341 Double_t AliKFParticleBase::GetSCorrection( const Double_t Part[], const Double_t XYZ[] ) 
342 {
343   //* Get big enough correction for S error to let the particle Part be fitted to XYZ point
344   
345   Double_t d[3] = { XYZ[0]-Part[0], XYZ[1]-Part[1], XYZ[2]-Part[2] };
346   Double_t p2 = Part[3]*Part[3]+Part[4]*Part[4]+Part[5]*Part[5];
347   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.;
348   return sigmaS;
349 }
350
351 void AliKFParticleBase::GetMeasurement( const Double_t XYZ[], Double_t m[], Double_t V[] ) const
352 {
353   //* Get additional covariances V used during measurement
354
355   Double_t b[3];
356   GetFieldValue( XYZ, b );
357   const Double_t kCLight =  0.000299792458;
358   b[0]*=kCLight; b[1]*=kCLight; b[2]*=kCLight;
359
360   Transport( GetDStoPoint(XYZ), m, V );
361
362   Double_t sigmaS = GetSCorrection( m, XYZ );
363
364   Double_t h[6];
365
366   h[0] = m[3]*sigmaS;
367   h[1] = m[4]*sigmaS;
368   h[2] = m[5]*sigmaS;
369   h[3] = ( h[1]*b[2]-h[2]*b[1] )*GetQ();
370   h[4] = ( h[2]*b[0]-h[0]*b[2] )*GetQ();
371   h[5] = ( h[0]*b[1]-h[1]*b[0] )*GetQ();
372     
373   V[ 0]+= h[0]*h[0];
374   V[ 1]+= h[1]*h[0];
375   V[ 2]+= h[1]*h[1];
376   V[ 3]+= h[2]*h[0];
377   V[ 4]+= h[2]*h[1];
378   V[ 5]+= h[2]*h[2];
379
380   V[ 6]+= h[3]*h[0];
381   V[ 7]+= h[3]*h[1];
382   V[ 8]+= h[3]*h[2];
383   V[ 9]+= h[3]*h[3];
384
385   V[10]+= h[4]*h[0];
386   V[11]+= h[4]*h[1];
387   V[12]+= h[4]*h[2];
388   V[13]+= h[4]*h[3];
389   V[14]+= h[4]*h[4];
390
391   V[15]+= h[5]*h[0];
392   V[16]+= h[5]*h[1];
393   V[17]+= h[5]*h[2];
394   V[18]+= h[5]*h[3];
395   V[19]+= h[5]*h[4];
396   V[20]+= h[5]*h[5];
397 }
398
399 void AliKFParticleBase::AddDaughter( const AliKFParticleBase &Daughter )
400 {
401   if( fNDF<-1 ){ // first daughter -> just copy
402     fNDF   = -1;
403     fQ     =  Daughter.GetQ();
404     for( Int_t i=0; i<7; i++) fP[i] = Daughter.fP[i];
405     for( Int_t i=0; i<28; i++) fC[i] = Daughter.fC[i];
406     fSFromDecay = 0;
407     fMassHypo = Daughter.fMassHypo;
408     SumDaughterMass = Daughter.SumDaughterMass;
409     return;
410   }
411
412   if(fConstructMethod == 0)
413     AddDaughterWithEnergyFit(Daughter);
414   else if(fConstructMethod == 1)
415     AddDaughterWithEnergyCalc(Daughter);
416   else if(fConstructMethod == 2)
417     AddDaughterWithEnergyFitMC(Daughter);
418
419   SumDaughterMass += Daughter.SumDaughterMass;
420   fMassHypo = -1;
421 }
422
423 void AliKFParticleBase::AddDaughterWithEnergyFit( const AliKFParticleBase &Daughter )
424 {
425   //* Energy considered as an independent veriable, fitted independently from momentum, without any constraints on mass
426
427   //* Add daughter 
428
429   TransportToDecayVertex();
430
431   Double_t b[3]; 
432   Int_t maxIter = 1;
433
434   if( !fIsLinearized ){
435     if( fNDF==-1 ){
436       Double_t ds, ds1;
437       GetDStoParticle(Daughter, ds, ds1);      
438       TransportToDS( ds );
439       Double_t m[8];
440       Double_t mCd[36];       
441       Daughter.Transport( ds1, m, mCd );    
442       fVtxGuess[0] = .5*( fP[0] + m[0] );
443       fVtxGuess[1] = .5*( fP[1] + m[1] );
444       fVtxGuess[2] = .5*( fP[2] + m[2] );
445     } else {
446       fVtxGuess[0] = fP[0];
447       fVtxGuess[1] = fP[1];
448       fVtxGuess[2] = fP[2]; 
449     }
450     maxIter = 3;
451   }
452
453   for( Int_t iter=0; iter<maxIter; iter++ ){
454
455     {
456       GetFieldValue( fVtxGuess, b );
457       const Double_t kCLight =  0.000299792458;
458       b[0]*=kCLight; b[1]*=kCLight; b[2]*=kCLight;
459     }
460
461     Double_t *ffP = fP, *ffC = fC, tmpP[8], tmpC[36];
462     if( fNDF==-1 ){            
463       GetMeasurement( fVtxGuess, tmpP, tmpC );
464       ffP = tmpP;
465       ffC = tmpC;
466     }
467
468     Double_t m[8], mV[36];
469
470     if( Daughter.fC[35]>0 ){
471       Daughter.GetMeasurement( fVtxGuess, m, mV );
472     } else {
473       for( Int_t i=0; i<8; i++ ) m[i] = Daughter.fP[i];
474       for( Int_t i=0; i<36; i++ ) mV[i] = Daughter.fC[i];
475     }
476     //*
477
478     Double_t mS[6];
479     {
480       Double_t mSi[6] = { ffC[0]+mV[0], 
481                           ffC[1]+mV[1], ffC[2]+mV[2], 
482                           ffC[3]+mV[3], ffC[4]+mV[4], ffC[5]+mV[5] };
483
484       mS[0] = mSi[2]*mSi[5] - mSi[4]*mSi[4];
485       mS[1] = mSi[3]*mSi[4] - mSi[1]*mSi[5];
486       mS[2] = mSi[0]*mSi[5] - mSi[3]*mSi[3];
487       mS[3] = mSi[1]*mSi[4] - mSi[2]*mSi[3];
488       mS[4] = mSi[1]*mSi[3] - mSi[0]*mSi[4];
489       mS[5] = mSi[0]*mSi[2] - mSi[1]*mSi[1];     
490
491       Double_t s = ( mSi[0]*mS[0] + mSi[1]*mS[1] + mSi[3]*mS[3] );      
492       s = ( TMath::Abs(s) > 1.E-20 )  ?1./s :0;   
493       mS[0]*=s;
494       mS[1]*=s;
495       mS[2]*=s;
496       mS[3]*=s;
497       mS[4]*=s;
498       mS[5]*=s;
499     }
500     //* Residual (measured - estimated)
501
502     Double_t zeta[3] = { m[0]-ffP[0], m[1]-ffP[1], m[2]-ffP[2] };    
503
504     //* CHt = CH' - D'
505
506     Double_t mCHt0[7], mCHt1[7], mCHt2[7];
507
508     mCHt0[0]=ffC[ 0] ;       mCHt1[0]=ffC[ 1] ;       mCHt2[0]=ffC[ 3] ;
509     mCHt0[1]=ffC[ 1] ;       mCHt1[1]=ffC[ 2] ;       mCHt2[1]=ffC[ 4] ;
510     mCHt0[2]=ffC[ 3] ;       mCHt1[2]=ffC[ 4] ;       mCHt2[2]=ffC[ 5] ;
511     mCHt0[3]=ffC[ 6]-mV[ 6]; mCHt1[3]=ffC[ 7]-mV[ 7]; mCHt2[3]=ffC[ 8]-mV[ 8];
512     mCHt0[4]=ffC[10]-mV[10]; mCHt1[4]=ffC[11]-mV[11]; mCHt2[4]=ffC[12]-mV[12];
513     mCHt0[5]=ffC[15]-mV[15]; mCHt1[5]=ffC[16]-mV[16]; mCHt2[5]=ffC[17]-mV[17];
514     mCHt0[6]=ffC[21]-mV[21]; mCHt1[6]=ffC[22]-mV[22]; mCHt2[6]=ffC[23]-mV[23];
515   
516     //* Kalman gain K = mCH'*S
517     
518     Double_t k0[7], k1[7], k2[7];
519     
520     for(Int_t i=0;i<7;++i){
521       k0[i] = mCHt0[i]*mS[0] + mCHt1[i]*mS[1] + mCHt2[i]*mS[3];
522       k1[i] = mCHt0[i]*mS[1] + mCHt1[i]*mS[2] + mCHt2[i]*mS[4];
523       k2[i] = mCHt0[i]*mS[3] + mCHt1[i]*mS[4] + mCHt2[i]*mS[5];
524     }
525
526    //* New estimation of the vertex position 
527
528     if( iter<maxIter-1 ){
529       for(Int_t i=0; i<3; ++i) 
530         fVtxGuess[i]= ffP[i] + k0[i]*zeta[0]+k1[i]*zeta[1]+k2[i]*zeta[2];
531       continue;
532     }
533
534     // last itearation -> update the particle
535
536     //* Add the daughter momentum to the particle momentum
537     
538     ffP[ 3] += m[ 3];
539     ffP[ 4] += m[ 4];
540     ffP[ 5] += m[ 5];
541     ffP[ 6] += m[ 6];
542   
543     ffC[ 9] += mV[ 9];
544     ffC[13] += mV[13];
545     ffC[14] += mV[14];
546     ffC[18] += mV[18];
547     ffC[19] += mV[19];
548     ffC[20] += mV[20];
549     ffC[24] += mV[24];
550     ffC[25] += mV[25];
551     ffC[26] += mV[26];
552     ffC[27] += mV[27];
553     
554  
555    //* New estimation of the vertex position r += K*zeta
556     
557     for(Int_t i=0;i<7;++i) 
558       fP[i] = ffP[i] + k0[i]*zeta[0] + k1[i]*zeta[1] + k2[i]*zeta[2];
559     
560     //* New covariance matrix C -= K*(mCH')'
561
562     for(Int_t i=0, k=0;i<7;++i){
563       for(Int_t j=0;j<=i;++j,++k){
564         fC[k] = ffC[k] - (k0[i]*mCHt0[j] + k1[i]*mCHt1[j] + k2[i]*mCHt2[j] );
565       }
566     }
567   
568     //* Calculate Chi^2 
569
570     fNDF  += 2;
571     fQ    +=  Daughter.GetQ();
572     fSFromDecay = 0;    
573     fChi2 += (mS[0]*zeta[0] + mS[1]*zeta[1] + mS[3]*zeta[2])*zeta[0]
574       +      (mS[1]*zeta[0] + mS[2]*zeta[1] + mS[4]*zeta[2])*zeta[1]
575       +      (mS[3]*zeta[0] + mS[4]*zeta[1] + mS[5]*zeta[2])*zeta[2];     
576
577   }
578 }
579
580 void AliKFParticleBase::AddDaughterWithEnergyCalc( const AliKFParticleBase &Daughter )
581 {
582   //* Energy considered as a dependent variable, calculated from the momentum and mass hypothesis
583
584   //* Add daughter 
585
586   TransportToDecayVertex();
587
588   Double_t b[3]; 
589   Int_t maxIter = 1;
590
591   if( !fIsLinearized ){
592     if( fNDF==-1 ){
593       Double_t ds, ds1;
594       GetDStoParticle(Daughter, ds, ds1);      
595       TransportToDS( ds );
596       Double_t m[8];
597       Double_t mCd[36];       
598       Daughter.Transport( ds1, m, mCd );    
599       fVtxGuess[0] = .5*( fP[0] + m[0] );
600       fVtxGuess[1] = .5*( fP[1] + m[1] );
601       fVtxGuess[2] = .5*( fP[2] + m[2] );
602     } else {
603       fVtxGuess[0] = fP[0];
604       fVtxGuess[1] = fP[1];
605       fVtxGuess[2] = fP[2]; 
606     }
607     maxIter = 3;
608   }
609
610   for( Int_t iter=0; iter<maxIter; iter++ ){
611
612     {
613       GetFieldValue( fVtxGuess, b );
614       const Double_t kCLight =  0.000299792458;
615       b[0]*=kCLight; b[1]*=kCLight; b[2]*=kCLight;
616     }
617
618     Double_t *ffP = fP, *ffC = fC, tmpP[8], tmpC[36];
619     if( fNDF==-1 ){            
620       GetMeasurement( fVtxGuess, tmpP, tmpC );
621       ffP = tmpP;
622       ffC = tmpC;
623     }
624
625     Double_t m[8], mV[36];
626
627     if( Daughter.fC[35]>0 ){
628       Daughter.GetMeasurement( fVtxGuess, m, mV );
629     } else {
630       for( Int_t i=0; i<8; i++ ) m[i] = Daughter.fP[i];
631       for( Int_t i=0; i<36; i++ ) mV[i] = Daughter.fC[i];
632     }
633
634     double Mass_mf2 = m[6]*m[6] - (m[3]*m[3] + m[4]*m[4] + m[5]*m[5]);
635     double Mass_rf2 = fP[6]*fP[6] - (fP[3]*fP[3] + fP[4]*fP[4] + fP[5]*fP[5]);
636
637     //*
638
639     Double_t mS[6];
640     {
641       Double_t mSi[6] = { ffC[0]+mV[0], 
642                           ffC[1]+mV[1], ffC[2]+mV[2], 
643                           ffC[3]+mV[3], ffC[4]+mV[4], ffC[5]+mV[5] };
644
645       mS[0] = mSi[2]*mSi[5] - mSi[4]*mSi[4];
646       mS[1] = mSi[3]*mSi[4] - mSi[1]*mSi[5];
647       mS[2] = mSi[0]*mSi[5] - mSi[3]*mSi[3];
648       mS[3] = mSi[1]*mSi[4] - mSi[2]*mSi[3];
649       mS[4] = mSi[1]*mSi[3] - mSi[0]*mSi[4];
650       mS[5] = mSi[0]*mSi[2] - mSi[1]*mSi[1];     
651
652       Double_t s = ( mSi[0]*mS[0] + mSi[1]*mS[1] + mSi[3]*mS[3] );      
653
654       s = ( s > 1.E-20 )  ?1./s :0;       
655       mS[0]*=s;
656       mS[1]*=s;
657       mS[2]*=s;
658       mS[3]*=s;
659       mS[4]*=s;
660       mS[5]*=s;
661     }
662
663     //* Residual (measured - estimated)
664
665     Double_t zeta[3] = { m[0]-ffP[0], m[1]-ffP[1], m[2]-ffP[2] };    
666
667     //* CHt = CH' - D'
668
669     Double_t mCHt0[6], mCHt1[6], mCHt2[6];
670
671     mCHt0[0]=ffC[ 0] ;       mCHt1[0]=ffC[ 1] ;       mCHt2[0]=ffC[ 3] ;
672     mCHt0[1]=ffC[ 1] ;       mCHt1[1]=ffC[ 2] ;       mCHt2[1]=ffC[ 4] ;
673     mCHt0[2]=ffC[ 3] ;       mCHt1[2]=ffC[ 4] ;       mCHt2[2]=ffC[ 5] ;
674     mCHt0[3]=ffC[ 6]-mV[ 6]; mCHt1[3]=ffC[ 7]-mV[ 7]; mCHt2[3]=ffC[ 8]-mV[ 8];
675     mCHt0[4]=ffC[10]-mV[10]; mCHt1[4]=ffC[11]-mV[11]; mCHt2[4]=ffC[12]-mV[12];
676     mCHt0[5]=ffC[15]-mV[15]; mCHt1[5]=ffC[16]-mV[16]; mCHt2[5]=ffC[17]-mV[17];
677
678     //* Kalman gain K = mCH'*S
679
680     Double_t k0[6], k1[6], k2[6];
681
682     for(Int_t i=0;i<6;++i){
683       k0[i] = mCHt0[i]*mS[0] + mCHt1[i]*mS[1] + mCHt2[i]*mS[3];
684       k1[i] = mCHt0[i]*mS[1] + mCHt1[i]*mS[2] + mCHt2[i]*mS[4];
685       k2[i] = mCHt0[i]*mS[3] + mCHt1[i]*mS[4] + mCHt2[i]*mS[5];
686     }
687
688    //* New estimation of the vertex position 
689
690     if( iter<maxIter-1 ){
691       for(Int_t i=0; i<3; ++i) 
692         fVtxGuess[i]= ffP[i] + k0[i]*zeta[0]+k1[i]*zeta[1]+k2[i]*zeta[2];
693       continue;
694     }
695
696    //* find mf and Vf - optimum value of the measurement and its covariance matrix
697     //* mVHt = V*H'
698     Double_t mVHt0[6], mVHt1[6], mVHt2[6];
699
700     mVHt0[0]= mV[ 0] ; mVHt1[0]= mV[ 1] ; mVHt2[0]= mV[ 3] ;
701     mVHt0[1]= mV[ 1] ; mVHt1[1]= mV[ 2] ; mVHt2[1]= mV[ 4] ;
702     mVHt0[2]= mV[ 3] ; mVHt1[2]= mV[ 4] ; mVHt2[2]= mV[ 5] ;
703     mVHt0[3]= mV[ 6] ; mVHt1[3]= mV[ 7] ; mVHt2[3]= mV[ 8] ;
704     mVHt0[4]= mV[10] ; mVHt1[4]= mV[11] ; mVHt2[4]= mV[12] ;
705     mVHt0[5]= mV[15] ; mVHt1[5]= mV[16] ; mVHt2[5]= mV[17] ;
706
707     //* Kalman gain Km = mCH'*S
708
709     Double_t km0[6], km1[6], km2[6];
710
711     for(Int_t i=0;i<6;++i){
712       km0[i] = mVHt0[i]*mS[0] + mVHt1[i]*mS[1] + mVHt2[i]*mS[3];
713       km1[i] = mVHt0[i]*mS[1] + mVHt1[i]*mS[2] + mVHt2[i]*mS[4];
714       km2[i] = mVHt0[i]*mS[3] + mVHt1[i]*mS[4] + mVHt2[i]*mS[5];
715     }
716
717     Double_t mf[7] = { m[0], m[1], m[2], m[3], m[4], m[5], m[6] };
718
719     for(Int_t i=0;i<6;++i) 
720       mf[i] = mf[i] - km0[i]*zeta[0] - km1[i]*zeta[1] - km2[i]*zeta[2];
721
722     Double_t E_mf = TMath::Sqrt( Mass_mf2 + (mf[3]*mf[3] + mf[4]*mf[4] + mf[5]*mf[5]) );
723
724     Double_t Vf[28];
725     for(Int_t iC=0; iC<28; iC++)
726       Vf[iC] = mV[iC];
727
728     //* hmf = d(E_mf)/d(mf)
729     Double_t hmf[7];
730     if( TMath::Abs(E_mf) < 1.e-10) hmf[3] = 0; else hmf[3] = mf[3]/E_mf;
731     if( TMath::Abs(E_mf) < 1.e-10) hmf[4] = 0; else hmf[4] = mf[4]/E_mf;
732     if( TMath::Abs(E_mf) < 1.e-10) hmf[5] = 0; else hmf[5] = mf[5]/E_mf;
733 //    if( TMath::Abs(E_mf) < 1.e-10) hmf[6] = 0; else hmf[6] = mf[6]/E_mf;
734     hmf[6] = 0;
735
736     for(Int_t i=0, k=0;i<6;++i){
737       for(Int_t j=0;j<=i;++j,++k){
738         Vf[k] = Vf[k] - (km0[i]*mVHt0[j] + km1[i]*mVHt1[j] + km2[i]*mVHt2[j] );
739       }
740     }
741     Double_t Vf24 = Vf[24], Vf25 = Vf[25], Vf26 = Vf[26];
742     Vf[21] = Vf[6 ]*hmf[3] + Vf[10]*hmf[4] + Vf[15]*hmf[5] + Vf[21]*hmf[6];
743     Vf[22] = Vf[7 ]*hmf[3] + Vf[11]*hmf[4] + Vf[16]*hmf[5] + Vf[22]*hmf[6];
744     Vf[23] = Vf[8 ]*hmf[3] + Vf[12]*hmf[4] + Vf[17]*hmf[5] + Vf[23]*hmf[6];
745     Vf[24] = Vf[9 ]*hmf[3] + Vf[13]*hmf[4] + Vf[18]*hmf[5] + Vf[24]*hmf[6];
746     Vf[25] = Vf[13]*hmf[3] + Vf[14]*hmf[4] + Vf[19]*hmf[5] + Vf[25]*hmf[6];
747     Vf[26] = Vf[18]*hmf[3] + Vf[19]*hmf[4] + Vf[20]*hmf[5] + Vf[26]*hmf[6];
748     Vf[27] = Vf[24]*hmf[3] + Vf[25]*hmf[4] + Vf[26]*hmf[5] + (Vf24*hmf[3] + Vf25*hmf[4] + Vf26*hmf[5] + Vf[27]*hmf[6])*hmf[6]; //here Vf[] are already modified
749
750     mf[6] = E_mf;
751
752     //* find rf and Cf - optimum value of the measurement and its covariance matrix
753
754     //* mCCHt = C*H'
755     Double_t mCCHt0[6], mCCHt1[6], mCCHt2[6];
756
757     mCCHt0[0]=ffC[ 0]; mCCHt1[0]=ffC[ 1]; mCCHt2[0]=ffC[ 3];
758     mCCHt0[1]=ffC[ 1]; mCCHt1[1]=ffC[ 2]; mCCHt2[1]=ffC[ 4];
759     mCCHt0[2]=ffC[ 3]; mCCHt1[2]=ffC[ 4]; mCCHt2[2]=ffC[ 5];
760     mCCHt0[3]=ffC[ 6]; mCCHt1[3]=ffC[ 7]; mCCHt2[3]=ffC[ 8];
761     mCCHt0[4]=ffC[10]; mCCHt1[4]=ffC[11]; mCCHt2[4]=ffC[12];
762     mCCHt0[5]=ffC[15]; mCCHt1[5]=ffC[16]; mCCHt2[5]=ffC[17];
763
764     //* Kalman gain Krf = mCH'*S
765
766     Double_t krf0[6], krf1[6], krf2[6];
767
768     for(Int_t i=0;i<6;++i){
769       krf0[i] = mCCHt0[i]*mS[0] + mCCHt1[i]*mS[1] + mCCHt2[i]*mS[3];
770       krf1[i] = mCCHt0[i]*mS[1] + mCCHt1[i]*mS[2] + mCCHt2[i]*mS[4];
771       krf2[i] = mCCHt0[i]*mS[3] + mCCHt1[i]*mS[4] + mCCHt2[i]*mS[5];
772     }
773     Double_t rf[7] = { ffP[0], ffP[1], ffP[2], ffP[3], ffP[4], ffP[5], ffP[6] };
774
775     for(Int_t i=0;i<6;++i) 
776       rf[i] = rf[i] + krf0[i]*zeta[0] + krf1[i]*zeta[1] + krf2[i]*zeta[2];
777
778     Double_t E_rf = TMath::Sqrt( Mass_rf2 + (rf[3]*rf[3] + rf[4]*rf[4] + rf[5]*rf[5]) );
779
780     Double_t Cf[28];
781     for(Int_t iC=0; iC<28; iC++)
782       Cf[iC] = ffC[iC];
783     //* hrf = d(Erf)/d(rf)
784     Double_t hrf[7];
785     if( TMath::Abs(E_rf) < 1.e-10) hrf[3] = 0; else hrf[3] = rf[3]/E_rf;
786     if( TMath::Abs(E_rf) < 1.e-10) hrf[4] = 0; else hrf[4] = rf[4]/E_rf;
787     if( TMath::Abs(E_rf) < 1.e-10) hrf[5] = 0; else hrf[5] = rf[5]/E_rf;
788 //    if( TMath::Abs(E_rf) < 1.e-10) hrf[6] = 0; else hrf[6] = rf[6]/E_rf;
789     hrf[6] = 0;
790
791     for(Int_t i=0, k=0;i<6;++i){
792       for(Int_t j=0;j<=i;++j,++k){
793         Cf[k] = Cf[k] - (krf0[i]*mCCHt0[j] + krf1[i]*mCCHt1[j] + krf2[i]*mCCHt2[j] );
794       }
795     }
796     Double_t Cf24 = Cf[24], Cf25 = Cf[25], Cf26 = Cf[26];
797     Cf[21] = Cf[6 ]*hrf[3] + Cf[10]*hrf[4] + Cf[15]*hrf[5] + Cf[21]*hrf[6];
798     Cf[22] = Cf[7 ]*hrf[3] + Cf[11]*hrf[4] + Cf[16]*hrf[5] + Cf[22]*hrf[6];
799     Cf[23] = Cf[8 ]*hrf[3] + Cf[12]*hrf[4] + Cf[17]*hrf[5] + Cf[23]*hrf[6];
800     Cf[24] = Cf[9 ]*hrf[3] + Cf[13]*hrf[4] + Cf[18]*hrf[5] + Cf[24]*hrf[6];
801     Cf[25] = Cf[13]*hrf[3] + Cf[14]*hrf[4] + Cf[19]*hrf[5] + Cf[25]*hrf[6];
802     Cf[26] = Cf[18]*hrf[3] + Cf[19]*hrf[4] + Cf[20]*hrf[5] + Cf[26]*hrf[6];
803     Cf[27] = Cf[24]*hrf[3] + Cf[25]*hrf[4] + Cf[26]*hrf[5] + (Cf24*hrf[3] + Cf25*hrf[4] + Cf26*hrf[5] + Cf[27]*hrf[6])*hrf[6]; //here Cf[] are already modified
804
805     for(Int_t iC=21; iC<28; iC++)
806     {
807       ffC[iC] = Cf[iC];
808       mV[iC]  = Vf[iC];
809     }
810
811     fP[6] = E_rf + E_mf;
812     rf[6] = E_rf;
813
814     //Double_t Dvv[3][3]; do not need this
815     Double_t Dvp[3][3];
816     Double_t Dpv[3][3];
817     Double_t Dpp[3][3];
818     Double_t De[7];
819
820     for(int i=0; i<3; i++)
821     {
822       for(int j=0; j<3; j++)
823       {
824         Dvp[i][j] = km0[i+3]*mCCHt0[j] + km1[i+3]*mCCHt1[j] + km2[i+3]*mCCHt2[j];
825         Dpv[i][j] = km0[i]*mCCHt0[j+3] + km1[i]*mCCHt1[j+3] + km2[i]*mCCHt2[j+3];
826         Dpp[i][j] = km0[i+3]*mCCHt0[j+3] + km1[i+3]*mCCHt1[j+3] + km2[i+3]*mCCHt2[j+3];
827       }
828     }
829
830     De[0] = hmf[3]*Dvp[0][0] + hmf[4]*Dvp[1][0] + hmf[5]*Dvp[2][0];
831     De[1] = hmf[3]*Dvp[0][1] + hmf[4]*Dvp[1][1] + hmf[5]*Dvp[2][1];
832     De[2] = hmf[3]*Dvp[0][2] + hmf[4]*Dvp[1][2] + hmf[5]*Dvp[2][2];
833     De[3] = hmf[3]*Dpp[0][0] + hmf[4]*Dpp[1][0] + hmf[5]*Dpp[2][0];
834     De[4] = hmf[3]*Dpp[0][1] + hmf[4]*Dpp[1][1] + hmf[5]*Dpp[2][1];
835     De[5] = hmf[3]*Dpp[0][2] + hmf[4]*Dpp[1][2] + hmf[5]*Dpp[2][2];
836     De[6] = 2*(De[3]*hrf[3] + De[4]*hrf[4] + De[5]*hrf[5]);
837
838     // last itearation -> update the particle
839
840     //* Add the daughter momentum to the particle momentum
841
842     ffP[ 3] += m[ 3];
843     ffP[ 4] += m[ 4];
844     ffP[ 5] += m[ 5];
845
846     ffC[ 9] += mV[ 9];
847     ffC[13] += mV[13];
848     ffC[14] += mV[14];
849     ffC[18] += mV[18];
850     ffC[19] += mV[19];
851     ffC[20] += mV[20];
852     ffC[24] += mV[24];
853     ffC[25] += mV[25];
854     ffC[26] += mV[26];
855     ffC[27] += mV[27];
856
857     ffC[21] += De[0];
858     ffC[22] += De[1];
859     ffC[23] += De[2];
860     ffC[24] += De[3];
861     ffC[25] += De[4];
862     ffC[26] += De[5];
863     ffC[27] += De[6];
864
865    //* New estimation of the vertex position r += K*zeta
866
867     for(Int_t i=0;i<6;++i) 
868       fP[i] = ffP[i] + k0[i]*zeta[0] + k1[i]*zeta[1] + k2[i]*zeta[2];
869
870     //* New covariance matrix C -= K*(mCH')'
871
872     for(Int_t i=0, k=0;i<6;++i){
873       for(Int_t j=0;j<=i;++j,++k){
874         fC[k] = ffC[k] - (k0[i]*mCHt0[j] + k1[i]*mCHt1[j] + k2[i]*mCHt2[j] );
875       }
876     }
877
878     for(int i=21; i<28; i++) fC[i] = ffC[i];
879
880     //* Calculate Chi^2 
881
882     fNDF  += 2;
883     fQ    +=  Daughter.GetQ();
884     fSFromDecay = 0;    
885     fChi2 += (mS[0]*zeta[0] + mS[1]*zeta[1] + mS[3]*zeta[2])*zeta[0]
886       +      (mS[1]*zeta[0] + mS[2]*zeta[1] + mS[4]*zeta[2])*zeta[1]
887       +      (mS[3]*zeta[0] + mS[4]*zeta[1] + mS[5]*zeta[2])*zeta[2];     
888   }
889 }
890
891 void AliKFParticleBase::AddDaughterWithEnergyFitMC( const AliKFParticleBase &Daughter )
892 {
893   //* Energy considered as an independent variable, fitted independently from momentum, without any constraints on mass
894
895   //* Add daughter 
896
897   TransportToDecayVertex();
898
899   Double_t b[3]; 
900   Int_t maxIter = 1;
901
902   if( !fIsLinearized ){
903     if( fNDF==-1 ){
904       Double_t ds, ds1;
905       GetDStoParticle(Daughter, ds, ds1);      
906       TransportToDS( ds );
907       Double_t m[8];
908       Double_t mCd[36];       
909       Daughter.Transport( ds1, m, mCd );    
910       fVtxGuess[0] = .5*( fP[0] + m[0] );
911       fVtxGuess[1] = .5*( fP[1] + m[1] );
912       fVtxGuess[2] = .5*( fP[2] + m[2] );
913     } else {
914       fVtxGuess[0] = fP[0];
915       fVtxGuess[1] = fP[1];
916       fVtxGuess[2] = fP[2]; 
917     }
918     maxIter = 3;
919   }
920
921   for( Int_t iter=0; iter<maxIter; iter++ ){
922
923     {
924       GetFieldValue( fVtxGuess, b );
925       const Double_t kCLight =  0.000299792458;
926       b[0]*=kCLight; b[1]*=kCLight; b[2]*=kCLight;
927     }
928
929     Double_t *ffP = fP, *ffC = fC, tmpP[8], tmpC[36];
930     if( fNDF==-1 ){            
931       GetMeasurement( fVtxGuess, tmpP, tmpC );
932       ffP = tmpP;
933       ffC = tmpC;
934     }
935     Double_t m[8], mV[36];
936
937     if( Daughter.fC[35]>0 ){
938       Daughter.GetMeasurement( fVtxGuess, m, mV );
939     } else {
940       for( Int_t i=0; i<8; i++ ) m[i] = Daughter.fP[i];
941       for( Int_t i=0; i<36; i++ ) mV[i] = Daughter.fC[i];
942     }
943     //*
944
945     Double_t mS[6];
946     {
947       Double_t mSi[6] = { ffC[0]+mV[0], 
948                           ffC[1]+mV[1], ffC[2]+mV[2], 
949                           ffC[3]+mV[3], ffC[4]+mV[4], ffC[5]+mV[5] };
950      
951       mS[0] = mSi[2]*mSi[5] - mSi[4]*mSi[4];
952       mS[1] = mSi[3]*mSi[4] - mSi[1]*mSi[5];
953       mS[2] = mSi[0]*mSi[5] - mSi[3]*mSi[3];
954       mS[3] = mSi[1]*mSi[4] - mSi[2]*mSi[3];
955       mS[4] = mSi[1]*mSi[3] - mSi[0]*mSi[4];
956       mS[5] = mSi[0]*mSi[2] - mSi[1]*mSi[1];     
957       
958       Double_t s = ( mSi[0]*mS[0] + mSi[1]*mS[1] + mSi[3]*mS[3] );      
959
960       s = ( s > 1.E-20 )  ?1./s :0;       
961       mS[0]*=s;
962       mS[1]*=s;
963       mS[2]*=s;
964       mS[3]*=s;
965       mS[4]*=s;
966       mS[5]*=s;
967     }
968     //* Residual (measured - estimated)
969     
970     Double_t zeta[3] = { m[0]-ffP[0], m[1]-ffP[1], m[2]-ffP[2] };    
971
972     
973     //* CHt = CH'
974     
975     Double_t mCHt0[7], mCHt1[7], mCHt2[7];
976     
977     mCHt0[0]=ffC[ 0] ; mCHt1[0]=ffC[ 1] ; mCHt2[0]=ffC[ 3] ;
978     mCHt0[1]=ffC[ 1] ; mCHt1[1]=ffC[ 2] ; mCHt2[1]=ffC[ 4] ;
979     mCHt0[2]=ffC[ 3] ; mCHt1[2]=ffC[ 4] ; mCHt2[2]=ffC[ 5] ;
980     mCHt0[3]=ffC[ 6] ; mCHt1[3]=ffC[ 7] ; mCHt2[3]=ffC[ 8] ;
981     mCHt0[4]=ffC[10] ; mCHt1[4]=ffC[11] ; mCHt2[4]=ffC[12] ;
982     mCHt0[5]=ffC[15] ; mCHt1[5]=ffC[16] ; mCHt2[5]=ffC[17] ;
983     mCHt0[6]=ffC[21] ; mCHt1[6]=ffC[22] ; mCHt2[6]=ffC[23] ;
984   
985     //* Kalman gain K = mCH'*S
986     
987     Double_t k0[7], k1[7], k2[7];
988     
989     for(Int_t i=0;i<7;++i){
990       k0[i] = mCHt0[i]*mS[0] + mCHt1[i]*mS[1] + mCHt2[i]*mS[3];
991       k1[i] = mCHt0[i]*mS[1] + mCHt1[i]*mS[2] + mCHt2[i]*mS[4];
992       k2[i] = mCHt0[i]*mS[3] + mCHt1[i]*mS[4] + mCHt2[i]*mS[5];
993     }
994
995    //* New estimation of the vertex position 
996
997     if( iter<maxIter-1 ){
998       for(Int_t i=0; i<3; ++i) 
999         fVtxGuess[i]= ffP[i] + k0[i]*zeta[0]+k1[i]*zeta[1]+k2[i]*zeta[2];
1000       continue;
1001     }
1002
1003     // last itearation -> update the particle
1004
1005     //* VHt = VH'
1006     
1007     Double_t mVHt0[7], mVHt1[7], mVHt2[7];
1008     
1009     mVHt0[0]=mV[ 0] ; mVHt1[0]=mV[ 1] ; mVHt2[0]=mV[ 3] ;
1010     mVHt0[1]=mV[ 1] ; mVHt1[1]=mV[ 2] ; mVHt2[1]=mV[ 4] ;
1011     mVHt0[2]=mV[ 3] ; mVHt1[2]=mV[ 4] ; mVHt2[2]=mV[ 5] ;
1012     mVHt0[3]=mV[ 6] ; mVHt1[3]=mV[ 7] ; mVHt2[3]=mV[ 8] ;
1013     mVHt0[4]=mV[10] ; mVHt1[4]=mV[11] ; mVHt2[4]=mV[12] ;
1014     mVHt0[5]=mV[15] ; mVHt1[5]=mV[16] ; mVHt2[5]=mV[17] ;
1015     mVHt0[6]=mV[21] ; mVHt1[6]=mV[22] ; mVHt2[6]=mV[23] ;
1016   
1017     //* Kalman gain Km = mCH'*S
1018     
1019     Double_t km0[7], km1[7], km2[7];
1020     
1021     for(Int_t i=0;i<7;++i){
1022       km0[i] = mVHt0[i]*mS[0] + mVHt1[i]*mS[1] + mVHt2[i]*mS[3];
1023       km1[i] = mVHt0[i]*mS[1] + mVHt1[i]*mS[2] + mVHt2[i]*mS[4];
1024       km2[i] = mVHt0[i]*mS[3] + mVHt1[i]*mS[4] + mVHt2[i]*mS[5];
1025     }
1026
1027     for(Int_t i=0;i<7;++i) 
1028       ffP[i] = ffP[i] + k0[i]*zeta[0] + k1[i]*zeta[1] + k2[i]*zeta[2];
1029
1030     for(Int_t i=0;i<7;++i) 
1031       m[i] = m[i] - km0[i]*zeta[0] - km1[i]*zeta[1] - km2[i]*zeta[2];
1032
1033     for(Int_t i=0, k=0;i<7;++i){
1034       for(Int_t j=0;j<=i;++j,++k){
1035         ffC[k] = ffC[k] - (k0[i]*mCHt0[j] + k1[i]*mCHt1[j] + k2[i]*mCHt2[j] );
1036       }
1037     }
1038
1039     for(Int_t i=0, k=0;i<7;++i){
1040       for(Int_t j=0;j<=i;++j,++k){
1041         mV[k] = mV[k] - (km0[i]*mVHt0[j] + km1[i]*mVHt1[j] + km2[i]*mVHt2[j] );
1042       }
1043     }
1044
1045     Double_t Df[7][7];
1046
1047     for(Int_t i=0;i<7;++i){
1048       for(Int_t j=0;j<7;++j){
1049         Df[i][j] = (km0[i]*mCHt0[j] + km1[i]*mCHt1[j] + km2[i]*mCHt2[j] );
1050       }
1051     }
1052
1053     Double_t J1[7][7], J2[7][7];
1054     for(Int_t iPar1=0; iPar1<7; iPar1++)
1055     {
1056       for(Int_t iPar2=0; iPar2<7; iPar2++)
1057       {
1058         J1[iPar1][iPar2] = 0;
1059         J2[iPar1][iPar2] = 0;
1060       }
1061     }
1062
1063     Double_t mMassParticle  = ffP[6]*ffP[6] - (ffP[3]*ffP[3] + ffP[4]*ffP[4] + ffP[5]*ffP[5]);
1064     Double_t mMassDaughter  = m[6]*m[6] - (m[3]*m[3] + m[4]*m[4] + m[5]*m[5]);
1065     if(mMassParticle > 0) mMassParticle = TMath::Sqrt(mMassParticle);
1066     if(mMassDaughter > 0) mMassDaughter = TMath::Sqrt(mMassDaughter);
1067
1068     if( fMassHypo > -0.5)
1069       SetMassConstraint(ffP,ffC,J1,fMassHypo);
1070     else if((mMassParticle < SumDaughterMass) || (ffP[6]<0) )
1071       SetMassConstraint(ffP,ffC,J1,SumDaughterMass);
1072
1073     if(Daughter.fMassHypo > -0.5)
1074       SetMassConstraint(m,mV,J2,Daughter.fMassHypo);
1075     else if((mMassDaughter < Daughter.SumDaughterMass) || (m[6] < 0) )
1076       SetMassConstraint(m,mV,J2,Daughter.SumDaughterMass);
1077
1078     Double_t DJ[7][7];
1079
1080     for(Int_t i=0; i<7; i++) {
1081       for(Int_t j=0; j<7; j++) {
1082         DJ[i][j] = 0;
1083         for(Int_t k=0; k<7; k++) {
1084           DJ[i][j] += Df[i][k]*J1[j][k];
1085         }
1086       }
1087     }
1088
1089     for(Int_t i=0; i<7; ++i){
1090       for(Int_t j=0; j<7; ++j){
1091         Df[i][j]=0;
1092         for(Int_t l=0; l<7; l++){
1093           Df[i][j] += J2[i][l]*DJ[l][j];
1094         }
1095       }
1096     }
1097
1098     //* Add the daughter momentum to the particle momentum
1099
1100     ffP[ 3] += m[ 3];
1101     ffP[ 4] += m[ 4];
1102     ffP[ 5] += m[ 5];
1103     ffP[ 6] += m[ 6];
1104
1105     ffC[ 9] += mV[ 9];
1106     ffC[13] += mV[13];
1107     ffC[14] += mV[14];
1108     ffC[18] += mV[18];
1109     ffC[19] += mV[19];
1110     ffC[20] += mV[20];
1111     ffC[24] += mV[24];
1112     ffC[25] += mV[25];
1113     ffC[26] += mV[26];
1114     ffC[27] += mV[27];
1115
1116     ffC[6 ] += Df[3][0]; ffC[7 ] += Df[3][1]; ffC[8 ] += Df[3][2];
1117     ffC[10] += Df[4][0]; ffC[11] += Df[4][1]; ffC[12] += Df[4][2];
1118     ffC[15] += Df[5][0]; ffC[16] += Df[5][1]; ffC[17] += Df[5][2];
1119     ffC[21] += Df[6][0]; ffC[22] += Df[6][1]; ffC[23] += Df[6][2];
1120
1121     ffC[9 ] += Df[3][3] + Df[3][3];
1122     ffC[13] += Df[4][3] + Df[3][4]; ffC[14] += Df[4][4] + Df[4][4];
1123     ffC[18] += Df[5][3] + Df[3][5]; ffC[19] += Df[5][4] + Df[4][5]; ffC[20] += Df[5][5] + Df[5][5];
1124     ffC[24] += Df[6][3] + Df[3][6]; ffC[25] += Df[6][4] + Df[4][6]; ffC[26] += Df[6][5] + Df[5][6]; ffC[27] += Df[6][6] + Df[6][6];
1125
1126    //* New estimation of the vertex position r += K*zeta
1127
1128     for(Int_t i=0;i<7;++i) 
1129       fP[i] = ffP[i];
1130
1131     //* New covariance matrix C -= K*(mCH')'
1132
1133     for(Int_t i=0, k=0;i<7;++i){
1134       for(Int_t j=0;j<=i;++j,++k){
1135         fC[k] = ffC[k];
1136       }
1137     }
1138     //* Calculate Chi^2 
1139
1140     fNDF  += 2;
1141     fQ    +=  Daughter.GetQ();
1142     fSFromDecay = 0;    
1143     fChi2 += (mS[0]*zeta[0] + mS[1]*zeta[1] + mS[3]*zeta[2])*zeta[0]
1144       +      (mS[1]*zeta[0] + mS[2]*zeta[1] + mS[4]*zeta[2])*zeta[1]
1145       +      (mS[3]*zeta[0] + mS[4]*zeta[1] + mS[5]*zeta[2])*zeta[2];
1146   }
1147 }
1148
1149 void AliKFParticleBase::SetProductionVertex( const AliKFParticleBase &Vtx )
1150 {
1151   //* Set production vertex for the particle, when the particle was not used in the vertex fit
1152
1153   const Double_t *m = Vtx.fP, *mV = Vtx.fC;
1154
1155   Bool_t noS = ( fC[35]<=0 ); // no decay length allowed
1156
1157   if( noS ){ 
1158     TransportToDecayVertex();
1159     fP[7] = 0;
1160     fC[28] = fC[29] = fC[30] = fC[31] = fC[32] = fC[33] = fC[34] = fC[35] = 0;
1161   } else {
1162     TransportToDS( GetDStoPoint( m ) );    
1163     fP[7] = -fSFromDecay;
1164     fC[28] = fC[29] = fC[30] = fC[31] = fC[32] = fC[33] = fC[34] = 0;
1165     fC[35] = 0.1;
1166     
1167     Convert(1);
1168   }
1169
1170   Double_t mAi[6];
1171
1172   InvertSym3( fC, mAi );
1173
1174   Double_t mB[5][3];
1175
1176   mB[0][0] = fC[ 6]*mAi[0] + fC[ 7]*mAi[1] + fC[ 8]*mAi[3];
1177   mB[0][1] = fC[ 6]*mAi[1] + fC[ 7]*mAi[2] + fC[ 8]*mAi[4];
1178   mB[0][2] = fC[ 6]*mAi[3] + fC[ 7]*mAi[4] + fC[ 8]*mAi[5];
1179
1180   mB[1][0] = fC[10]*mAi[0] + fC[11]*mAi[1] + fC[12]*mAi[3];
1181   mB[1][1] = fC[10]*mAi[1] + fC[11]*mAi[2] + fC[12]*mAi[4];
1182   mB[1][2] = fC[10]*mAi[3] + fC[11]*mAi[4] + fC[12]*mAi[5];
1183
1184   mB[2][0] = fC[15]*mAi[0] + fC[16]*mAi[1] + fC[17]*mAi[3];
1185   mB[2][1] = fC[15]*mAi[1] + fC[16]*mAi[2] + fC[17]*mAi[4];
1186   mB[2][2] = fC[15]*mAi[3] + fC[16]*mAi[4] + fC[17]*mAi[5];
1187
1188   mB[3][0] = fC[21]*mAi[0] + fC[22]*mAi[1] + fC[23]*mAi[3];
1189   mB[3][1] = fC[21]*mAi[1] + fC[22]*mAi[2] + fC[23]*mAi[4];
1190   mB[3][2] = fC[21]*mAi[3] + fC[22]*mAi[4] + fC[23]*mAi[5];
1191
1192   mB[4][0] = fC[28]*mAi[0] + fC[29]*mAi[1] + fC[30]*mAi[3];
1193   mB[4][1] = fC[28]*mAi[1] + fC[29]*mAi[2] + fC[30]*mAi[4];
1194   mB[4][2] = fC[28]*mAi[3] + fC[29]*mAi[4] + fC[30]*mAi[5];
1195
1196   Double_t z[3] = { m[0]-fP[0], m[1]-fP[1], m[2]-fP[2] };
1197
1198   {
1199     Double_t mAVi[6] = { fC[0]-mV[0], fC[1]-mV[1], fC[2]-mV[2], 
1200                         fC[3]-mV[3], fC[4]-mV[4], fC[5]-mV[5] };
1201     
1202     if( !InvertSym3( mAVi, mAVi ) ){
1203
1204       Double_t dChi2 = ( +(mAVi[0]*z[0] + mAVi[1]*z[1] + mAVi[3]*z[2])*z[0]
1205                          +(mAVi[1]*z[0] + mAVi[2]*z[1] + mAVi[4]*z[2])*z[1]
1206                          +(mAVi[3]*z[0] + mAVi[4]*z[1] + mAVi[5]*z[2])*z[2] );
1207       
1208       // Take Abs(dChi2) here. Negative value of 'det' or 'dChi2' shows that the particle 
1209       // was not used in the production vertex fit
1210       
1211       fChi2+= TMath::Abs( dChi2 );
1212     }
1213     fNDF  += 2;
1214   }
1215   
1216   fP[0] = m[0];
1217   fP[1] = m[1];
1218   fP[2] = m[2];
1219   fP[3]+= mB[0][0]*z[0] + mB[0][1]*z[1] + mB[0][2]*z[2];
1220   fP[4]+= mB[1][0]*z[0] + mB[1][1]*z[1] + mB[1][2]*z[2];
1221   fP[5]+= mB[2][0]*z[0] + mB[2][1]*z[1] + mB[2][2]*z[2];
1222   fP[6]+= mB[3][0]*z[0] + mB[3][1]*z[1] + mB[3][2]*z[2];
1223   fP[7]+= mB[4][0]*z[0] + mB[4][1]*z[1] + mB[4][2]*z[2];
1224   
1225   Double_t d0, d1, d2;
1226
1227   fC[0] = mV[0];
1228   fC[1] = mV[1];
1229   fC[2] = mV[2];
1230   fC[3] = mV[3];
1231   fC[4] = mV[4];
1232   fC[5] = mV[5];
1233
1234   d0= mB[0][0]*mV[0] + mB[0][1]*mV[1] + mB[0][2]*mV[3] - fC[ 6];
1235   d1= mB[0][0]*mV[1] + mB[0][1]*mV[2] + mB[0][2]*mV[4] - fC[ 7];
1236   d2= mB[0][0]*mV[3] + mB[0][1]*mV[4] + mB[0][2]*mV[5] - fC[ 8];
1237
1238   fC[ 6]+= d0;
1239   fC[ 7]+= d1;
1240   fC[ 8]+= d2;
1241   fC[ 9]+= d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
1242
1243   d0= mB[1][0]*mV[0] + mB[1][1]*mV[1] + mB[1][2]*mV[3] - fC[10];
1244   d1= mB[1][0]*mV[1] + mB[1][1]*mV[2] + mB[1][2]*mV[4] - fC[11];
1245   d2= mB[1][0]*mV[3] + mB[1][1]*mV[4] + mB[1][2]*mV[5] - fC[12];
1246
1247   fC[10]+= d0;
1248   fC[11]+= d1;
1249   fC[12]+= d2;
1250   fC[13]+= d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
1251   fC[14]+= d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
1252
1253   d0= mB[2][0]*mV[0] + mB[2][1]*mV[1] + mB[2][2]*mV[3] - fC[15];
1254   d1= mB[2][0]*mV[1] + mB[2][1]*mV[2] + mB[2][2]*mV[4] - fC[16];
1255   d2= mB[2][0]*mV[3] + mB[2][1]*mV[4] + mB[2][2]*mV[5] - fC[17];
1256
1257   fC[15]+= d0;
1258   fC[16]+= d1;
1259   fC[17]+= d2;
1260   fC[18]+= d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
1261   fC[19]+= d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
1262   fC[20]+= d0*mB[2][0] + d1*mB[2][1] + d2*mB[2][2];
1263
1264   d0= mB[3][0]*mV[0] + mB[3][1]*mV[1] + mB[3][2]*mV[3] - fC[21];
1265   d1= mB[3][0]*mV[1] + mB[3][1]*mV[2] + mB[3][2]*mV[4] - fC[22];
1266   d2= mB[3][0]*mV[3] + mB[3][1]*mV[4] + mB[3][2]*mV[5] - fC[23];
1267
1268   fC[21]+= d0;
1269   fC[22]+= d1;
1270   fC[23]+= d2;
1271   fC[24]+= d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
1272   fC[25]+= d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
1273   fC[26]+= d0*mB[2][0] + d1*mB[2][1] + d2*mB[2][2];
1274   fC[27]+= d0*mB[3][0] + d1*mB[3][1] + d2*mB[3][2];
1275
1276   d0= mB[4][0]*mV[0] + mB[4][1]*mV[1] + mB[4][2]*mV[3] - fC[28];
1277   d1= mB[4][0]*mV[1] + mB[4][1]*mV[2] + mB[4][2]*mV[4] - fC[29];
1278   d2= mB[4][0]*mV[3] + mB[4][1]*mV[4] + mB[4][2]*mV[5] - fC[30];
1279
1280   fC[28]+= d0;
1281   fC[29]+= d1;
1282   fC[30]+= d2;
1283   fC[31]+= d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
1284   fC[32]+= d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
1285   fC[33]+= d0*mB[2][0] + d1*mB[2][1] + d2*mB[2][2];
1286   fC[34]+= d0*mB[3][0] + d1*mB[3][1] + d2*mB[3][2];
1287   fC[35]+= d0*mB[4][0] + d1*mB[4][1] + d2*mB[4][2];
1288   
1289   if( noS ){ 
1290     fP[7] = 0;
1291     fC[28] = fC[29] = fC[30] = fC[31] = fC[32] = fC[33] = fC[34] = fC[35] = 0;
1292   } else {
1293     TransportToDS( fP[7] );
1294     Convert(0);
1295   }
1296
1297   fSFromDecay = 0;
1298 }
1299
1300 void AliKFParticleBase::SetMassConstraint( Double_t *mP, Double_t *mC, Double_t J[7][7], Double_t Mass )
1301 {
1302   const Double_t E2 = mP[6]*mP[6], p2 = mP[3]*mP[3]+mP[4]*mP[4]+mP[5]*mP[5], M2 = Mass*Mass;
1303
1304   const Double_t a = E2 - p2 + 2.*M2;
1305   const Double_t b = -2.*(E2 + p2);
1306   const Double_t c = E2 - p2 - M2;
1307
1308   Double_t Lambda = 0;
1309   if(TMath::Abs(b) > 1.e-10) Lambda = -c / b;
1310
1311   Double_t D = 4.*E2*p2 - M2*(E2-p2-2.*M2);
1312   if(D>=0 && TMath::Abs(a) > 1.e-10) Lambda = (E2 + p2 - sqrt(D))/a;
1313
1314   if(mP[6] < 0) //If energy < 0 we need a Lambda < 0
1315     Lambda = -1000000.; //Empirical, a better solution should be found
1316
1317   Int_t iIter=0;
1318   for(iIter=0; iIter<100; iIter++)
1319   {
1320     Double_t Lambda2 = Lambda*Lambda;
1321     Double_t Lambda4 = Lambda2*Lambda2;
1322
1323     Double_t Lambda0 = Lambda;
1324
1325     Double_t f  = -M2 * Lambda4 + a*Lambda2 + b*Lambda + c;
1326     Double_t df = -4.*M2 * Lambda2*Lambda + 2.*a*Lambda + b;
1327     if(TMath::Abs(df) < 1.e-10) break;
1328     Lambda -= f/df;
1329     if(TMath::Abs(Lambda0 - Lambda) < 1.e-8) break;
1330   }
1331
1332   const Double_t LPi = 1./(1. + Lambda);
1333   const Double_t LMi = 1./(1. - Lambda);
1334   const Double_t LP2i = LPi*LPi;
1335   const Double_t LM2i = LMi*LMi;
1336
1337   Double_t Lambda2 = Lambda*Lambda;
1338
1339   Double_t dfl  = -4.*M2 * Lambda2*Lambda + 2.*a*Lambda + b;
1340   Double_t dfx[4] = {0,0,0,0};
1341   dfx[0] = -2.*(1. + Lambda)*(1. + Lambda)*mP[3];
1342   dfx[1] = -2.*(1. + Lambda)*(1. + Lambda)*mP[4];
1343   dfx[2] = -2.*(1. + Lambda)*(1. + Lambda)*mP[5];
1344   dfx[3] = 2.*(1. - Lambda)*(1. - Lambda)*mP[6];
1345   Double_t dlx[7] = {1,1,1,1,1,1,1};
1346   if(TMath::Abs(dfl) > 1.e-10 )
1347   {
1348     for(int i=0; i<7; i++)
1349       dlx[i] = -dfx[i] / dfl;
1350   }
1351
1352   Double_t dxx[4] = {mP[3]*LM2i, mP[4]*LM2i, mP[5]*LM2i, -mP[6]*LP2i};
1353
1354   for(Int_t i=0; i<7; i++)
1355     for(Int_t j=0; j<7; j++)
1356       J[i][j]=0;
1357   J[0][0] = 1.;
1358   J[1][1] = 1.;
1359   J[2][2] = 1.;
1360
1361   for(Int_t i=3; i<7; i++)
1362     for(Int_t j=3; j<7; j++)
1363       J[i][j] = dlx[j-3]*dxx[i-3];
1364
1365   for(Int_t i=3; i<6; i++)
1366     J[i][i] += LMi;
1367   J[6][6] += LPi;
1368
1369   Double_t CJ[7][7];
1370
1371   for(Int_t i=0; i<7; i++) {
1372     for(Int_t j=0; j<7; j++) {
1373       CJ[i][j] = 0;
1374       for(Int_t k=0; k<7; k++) {
1375         CJ[i][j] += mC[IJ(i,k)]*J[j][k];
1376       }
1377     }
1378   }
1379
1380   for(Int_t i=0; i<7; ++i){
1381     for(Int_t j=0; j<=i; ++j){
1382       mC[IJ(i,j)]=0;
1383       for(Int_t l=0; l<7; l++){
1384         mC[IJ(i,j)] += J[i][l]*CJ[l][j];
1385       }
1386     }
1387   }
1388
1389   mP[3] *= LMi;
1390   mP[4] *= LMi;
1391   mP[5] *= LMi;
1392   mP[6] *= LPi;
1393 }
1394
1395 void AliKFParticleBase::SetNonlinearMassConstraint( Double_t Mass )
1396 {
1397   Double_t J[7][7];
1398   SetMassConstraint( fP, fC, J, Mass );
1399   fMassHypo = Mass;
1400   SumDaughterMass = Mass;
1401 }
1402
1403 void AliKFParticleBase::SetMassConstraint( Double_t Mass, Double_t SigmaMass )
1404 {  
1405   //* Set hard( SigmaMass=0 ) or soft (SigmaMass>0) mass constraint 
1406
1407   fMassHypo = Mass;
1408   SumDaughterMass = Mass;
1409
1410   Double_t m2 = Mass*Mass;            // measurement, weighted by Mass 
1411   Double_t s2 = m2*SigmaMass*SigmaMass; // sigma^2
1412
1413   Double_t p2 = fP[3]*fP[3] + fP[4]*fP[4] + fP[5]*fP[5]; 
1414   Double_t e0 = TMath::Sqrt(m2+p2);
1415
1416   Double_t mH[8];
1417   mH[0] = mH[1] = mH[2] = 0.;
1418   mH[3] = -2*fP[3]; 
1419   mH[4] = -2*fP[4]; 
1420   mH[5] = -2*fP[5]; 
1421   mH[6] =  2*fP[6];//e0;
1422   mH[7] = 0; 
1423
1424   Double_t zeta = e0*e0 - e0*fP[6];
1425   zeta = m2 - (fP[6]*fP[6]-p2);
1426   
1427   Double_t mCHt[8], s2_est=0;
1428   for( Int_t i=0; i<8; ++i ){
1429     mCHt[i] = 0.0;
1430     for (Int_t j=0;j<8;++j) mCHt[i] += Cij(i,j)*mH[j];
1431     s2_est += mH[i]*mCHt[i];
1432   }
1433   
1434   if( s2_est<1.e-20 ) return; // calculated mass error is already 0, 
1435                               // the particle can not be constrained on mass
1436
1437   Double_t w2 = 1./( s2 + s2_est );
1438   fChi2 += zeta*zeta*w2;
1439   fNDF  += 1;
1440   for( Int_t i=0, ii=0; i<8; ++i ){
1441     Double_t ki = mCHt[i]*w2;
1442     fP[i]+= ki*zeta;
1443     for(Int_t j=0;j<=i;++j) fC[ii++] -= ki*mCHt[j];    
1444   }
1445 }
1446
1447
1448 void AliKFParticleBase::SetNoDecayLength()
1449 {  
1450   //* Set no decay length for resonances
1451
1452   TransportToDecayVertex();
1453
1454   Double_t h[8];
1455   h[0] = h[1] = h[2] = h[3] = h[4] = h[5] = h[6] = 0;
1456   h[7] = 1; 
1457
1458   Double_t zeta = 0 - fP[7];
1459   for(Int_t i=0;i<8;++i) zeta -= h[i]*(fP[i]-fP[i]);
1460   
1461   Double_t s = fC[35];   
1462   if( s>1.e-20 ){
1463     s = 1./s;
1464     fChi2 += zeta*zeta*s;
1465     fNDF  += 1;
1466     for( Int_t i=0, ii=0; i<7; ++i ){
1467       Double_t ki = fC[28+i]*s;
1468       fP[i]+= ki*zeta;
1469       for(Int_t j=0;j<=i;++j) fC[ii++] -= ki*fC[28+j];    
1470     }
1471   }
1472   fP[7] = 0;
1473   fC[28] = fC[29] = fC[30] = fC[31] = fC[32] = fC[33] = fC[34] = fC[35] = 0;
1474 }
1475
1476
1477 void AliKFParticleBase::Construct( const AliKFParticleBase* vDaughters[], Int_t NDaughters,
1478                                    const AliKFParticleBase *Parent,  Double_t Mass, Bool_t IsConstrained         )
1479
1480   //* Full reconstruction in one go
1481
1482   Int_t maxIter = 1;
1483   bool wasLinearized = fIsLinearized;
1484   if( !fIsLinearized || IsConstrained ){
1485     //fVtxGuess[0] = fVtxGuess[1] = fVtxGuess[2] = 0;  //!!!!
1486     fVtxGuess[0] = GetX();
1487     fVtxGuess[1] = GetY();
1488     fVtxGuess[2] = GetZ();
1489     fIsLinearized = 1;
1490     maxIter = 3;
1491   }
1492
1493   Double_t constraintC[6];
1494
1495   if( IsConstrained ){
1496     for(Int_t i=0;i<6;++i) constraintC[i]=fC[i];
1497   } else {
1498     for(Int_t i=0;i<6;++i) constraintC[i]=0.;
1499     constraintC[0] = constraintC[2] = constraintC[5] = 100.;    
1500   }
1501
1502
1503   for( Int_t iter=0; iter<maxIter; iter++ ){
1504     fAtProductionVertex = 0;
1505     fSFromDecay = 0;
1506     fP[0] = fVtxGuess[0];
1507     fP[1] = fVtxGuess[1];
1508     fP[2] = fVtxGuess[2];
1509     fP[3] = 0;
1510     fP[4] = 0;
1511     fP[5] = 0;
1512     fP[6] = 0;
1513     fP[7] = 0;
1514     SumDaughterMass = 0;
1515
1516     for(Int_t i=0;i<6; ++i) fC[i]=constraintC[i];
1517     for(Int_t i=6;i<36;++i) fC[i]=0.;
1518     fC[35] = 1.;
1519
1520     fNDF  = IsConstrained ?0 :-3;
1521     fChi2 =  0.;
1522     fQ = 0;
1523
1524     for( Int_t itr =0; itr<NDaughters; itr++ ){
1525       AddDaughter( *vDaughters[itr] );    
1526     }
1527     if( iter<maxIter-1){
1528       for( Int_t i=0; i<3; i++ ) fVtxGuess[i] = fP[i];  
1529     }
1530   }
1531   fIsLinearized = wasLinearized;    
1532
1533   if( Mass>=0 ) SetMassConstraint( Mass );
1534   if( Parent ) SetProductionVertex( *Parent );
1535 }
1536
1537
1538 void AliKFParticleBase::Convert( bool ToProduction )
1539 {
1540   //* Tricky function - convert the particle error along its trajectory to 
1541   //* the value which corresponds to its production/decay vertex
1542   //* It is done by combination of the error of decay length with the position errors
1543
1544   Double_t fld[3]; 
1545   {
1546     GetFieldValue( fP, fld );
1547     const Double_t kCLight =  fQ*0.000299792458;
1548     fld[0]*=kCLight; fld[1]*=kCLight; fld[2]*=kCLight;
1549   }
1550
1551   Double_t h[6];
1552   
1553   h[0] = fP[3];
1554   h[1] = fP[4];
1555   h[2] = fP[5];
1556   if( ToProduction ){ h[0]=-h[0]; h[1]=-h[1]; h[2]=-h[2]; } 
1557   h[3] = h[1]*fld[2]-h[2]*fld[1];
1558   h[4] = h[2]*fld[0]-h[0]*fld[2];
1559   h[5] = h[0]*fld[1]-h[1]*fld[0];
1560   
1561   Double_t c;
1562
1563   c = fC[28]+h[0]*fC[35];
1564   fC[ 0]+= h[0]*(c+fC[28]);
1565   fC[28] = c;
1566
1567   fC[ 1]+= h[1]*fC[28] + h[0]*fC[29];
1568   c = fC[29]+h[1]*fC[35];
1569   fC[ 2]+= h[1]*(c+fC[29]);
1570   fC[29] = c;
1571
1572   fC[ 3]+= h[2]*fC[28] + h[0]*fC[30];
1573   fC[ 4]+= h[2]*fC[29] + h[1]*fC[30];
1574   c = fC[30]+h[2]*fC[35];
1575   fC[ 5]+= h[2]*(c+fC[30]);
1576   fC[30] = c;
1577
1578   fC[ 6]+= h[3]*fC[28] + h[0]*fC[31];
1579   fC[ 7]+= h[3]*fC[29] + h[1]*fC[31];
1580   fC[ 8]+= h[3]*fC[30] + h[2]*fC[31];
1581   c = fC[31]+h[3]*fC[35];
1582   fC[ 9]+= h[3]*(c+fC[31]);
1583   fC[31] = c;
1584   
1585   fC[10]+= h[4]*fC[28] + h[0]*fC[32];
1586   fC[11]+= h[4]*fC[29] + h[1]*fC[32];
1587   fC[12]+= h[4]*fC[30] + h[2]*fC[32];
1588   fC[13]+= h[4]*fC[31] + h[3]*fC[32];
1589   c = fC[32]+h[4]*fC[35];
1590   fC[14]+= h[4]*(c+fC[32]);
1591   fC[32] = c;
1592   
1593   fC[15]+= h[5]*fC[28] + h[0]*fC[33];
1594   fC[16]+= h[5]*fC[29] + h[1]*fC[33];
1595   fC[17]+= h[5]*fC[30] + h[2]*fC[33];
1596   fC[18]+= h[5]*fC[31] + h[3]*fC[33];
1597   fC[19]+= h[5]*fC[32] + h[4]*fC[33];
1598   c = fC[33]+h[5]*fC[35];
1599   fC[20]+= h[5]*(c+fC[33]);
1600   fC[33] = c;
1601
1602   fC[21]+= h[0]*fC[34];
1603   fC[22]+= h[1]*fC[34];
1604   fC[23]+= h[2]*fC[34];
1605   fC[24]+= h[3]*fC[34];
1606   fC[25]+= h[4]*fC[34];
1607   fC[26]+= h[5]*fC[34];
1608 }
1609
1610
1611 void AliKFParticleBase::TransportToDecayVertex()
1612 {
1613   //* Transport the particle to its decay vertex 
1614
1615   if( fSFromDecay != 0 ) TransportToDS( -fSFromDecay );
1616   if( fAtProductionVertex ) Convert(0);
1617   fAtProductionVertex = 0;
1618 }
1619
1620 void AliKFParticleBase::TransportToProductionVertex()
1621 {
1622   //* Transport the particle to its production vertex 
1623   
1624   if( fSFromDecay != -fP[7] ) TransportToDS( -fSFromDecay-fP[7] );
1625   if( !fAtProductionVertex ) Convert( 1 );
1626   fAtProductionVertex = 1;
1627 }
1628
1629
1630 void AliKFParticleBase::TransportToDS( Double_t dS )
1631
1632   //* Transport the particle on dS parameter (SignedPath/Momentum) 
1633  
1634   Transport( dS, fP, fC );
1635   fSFromDecay+= dS;
1636 }
1637
1638
1639 Double_t AliKFParticleBase::GetDStoPointLine( const Double_t xyz[] ) const 
1640 {
1641   //* Get dS to a certain space point without field
1642
1643   Double_t p2 = fP[3]*fP[3] + fP[4]*fP[4] + fP[5]*fP[5];  
1644   if( p2<1.e-4 ) p2 = 1;
1645   return ( fP[3]*(xyz[0]-fP[0]) + fP[4]*(xyz[1]-fP[1]) + fP[5]*(xyz[2]-fP[2]) )/p2;
1646 }
1647
1648
1649 Double_t AliKFParticleBase::GetDStoPointBz( Double_t B, const Double_t xyz[] ) 
1650   const
1651
1652   
1653   //* Get dS to a certain space point for Bz field
1654   const Double_t kCLight = 0.000299792458;
1655   Double_t bq = B*fQ*kCLight;
1656   Double_t pt2 = fP[3]*fP[3] + fP[4]*fP[4];
1657   if( pt2<1.e-4 ) return 0;
1658   Double_t dx = xyz[0] - fP[0];
1659   Double_t dy = xyz[1] - fP[1]; 
1660   Double_t a = dx*fP[3]+dy*fP[4];
1661   Double_t dS;
1662
1663   if( TMath::Abs(bq)<1.e-8 ) dS = a/pt2;  
1664   else dS =  TMath::ATan2( bq*a, pt2 + bq*(dy*fP[3] -dx*fP[4]) )/bq;
1665
1666   if(0){
1667
1668     Double_t px = fP[3];
1669     Double_t py = fP[4];
1670     Double_t pz = fP[5];
1671     Double_t ss[2], g[2][5];
1672   
1673     ss[0] = dS;
1674     ss[1] = -dS;
1675     for( Int_t i=0; i<2; i++){
1676       Double_t bs = bq*ss[i];
1677       Double_t c = TMath::Cos(bs), s = TMath::Sin(bs);
1678       Double_t cB,sB;
1679       if( TMath::Abs(bq)>1.e-8){
1680         cB= (1-c)/bq;     
1681         sB= s/bq;  
1682       }else{
1683         const Double_t kOvSqr6 = 1./TMath::Sqrt(6.);
1684         sB = (1.-bs*kOvSqr6)*(1.+bs*kOvSqr6)*ss[i];
1685         cB = .5*sB*bs;
1686       }
1687       g[i][0] = fP[0] + sB*px + cB*py;
1688       g[i][1] = fP[1] - cB*px + sB*py;
1689       g[i][2] = fP[2] + ss[i]*pz;
1690       g[i][3] =       + c*px + s*py;
1691       g[i][4] =       - s*px + c*py;      
1692     }
1693
1694     Int_t i=0;
1695   
1696     Double_t dMin = 1.e10;
1697     for( Int_t j=0; j<2; j++){
1698       Double_t xx = g[j][0]-xyz[0];
1699       Double_t yy = g[j][1]-xyz[1];
1700       Double_t zz = g[j][2]-xyz[2];
1701       Double_t d = xx*xx + yy*yy + zz*zz;
1702       if( d<dMin ){
1703         dMin = d;
1704         i = j;
1705       }
1706     }
1707
1708     dS = ss[i];
1709
1710     Double_t x= g[i][0], y= g[i][1], z= g[i][2], ppx= g[i][3], ppy= g[i][4];      
1711     Double_t ddx = x-xyz[0];
1712     Double_t ddy = y-xyz[1];
1713     Double_t ddz = z-xyz[2];
1714     Double_t c = ddx*ppx  + ddy*ppy  + ddz*pz ;
1715     Double_t pp2 = ppx*ppx + ppy*ppy + pz*pz;    
1716     if( TMath::Abs(pp2)>1.e-8 ){
1717       dS+=c/pp2;
1718     }
1719   }
1720   return dS;
1721 }
1722
1723
1724 void AliKFParticleBase::GetDStoParticleBz( Double_t B, const AliKFParticleBase &p, 
1725                                            Double_t &DS, Double_t &DS1 ) 
1726   const
1727
1728   //* Get dS to another particle for Bz field
1729   Double_t px = fP[3];
1730   Double_t py = fP[4];
1731   Double_t pz = fP[5];
1732
1733   Double_t px1 = p.fP[3];
1734   Double_t py1 = p.fP[4];
1735   Double_t pz1 = p.fP[5];
1736
1737   const Double_t kCLight = 0.000299792458;
1738
1739   Double_t bq = B*fQ*kCLight;
1740   Double_t bq1 = B*p.fQ*kCLight;
1741   Double_t s=0, ds=0, s1=0, ds1=0;
1742   
1743   if( TMath::Abs(bq)>1.e-8 || TMath::Abs(bq1)>1.e-8 ){
1744
1745     Double_t dx = (p.fP[0] - fP[0]);
1746     Double_t dy = (p.fP[1] - fP[1]);
1747     Double_t d2 = (dx*dx+dy*dy);
1748     
1749     Double_t p2  = (px *px  + py *py); 
1750     Double_t p21 = (px1*px1 + py1*py1);
1751
1752     if( TMath::Abs(p2) < 1.e-8 || TMath::Abs(p21) < 1.e-8 )
1753     {
1754       DS=0.;
1755       DS1=0.;
1756       return;
1757     }
1758
1759     Double_t a = (px*py1 - py*px1);
1760     Double_t b = (px*px1 + py*py1);
1761     
1762     Double_t ldx = bq*bq1*dx - bq1*py + bq*py1 ;
1763     Double_t ldy = bq*bq1*dy + bq1*px - bq*px1 ;
1764     Double_t l2 = ldx*ldx + ldy*ldy;
1765     
1766     Double_t cS = bq1*p2 + bq*bq1*(dy* px - dx* py) -  bq*b;
1767     Double_t cS1= bq*p21 - bq*bq1*(dy*px1 - dx*py1) - bq1*b;
1768
1769     Double_t ca  = bq*bq*bq1*d2  +2*( cS + bq*bq*(py1*dx-px1*dy)) ;
1770     Double_t ca1 = bq*bq1*bq1*d2 +2*( cS1 - bq1*bq1*(py*dx-px*dy)) ;  
1771   
1772     Double_t sa = 4*l2*p2 - ca*ca;
1773     Double_t sa1 = 4*l2*p21 - ca1*ca1;
1774
1775     if(sa<0) sa=0;
1776     if(sa1<0)sa1=0;
1777
1778     if( TMath::Abs(bq)>1.e-8){
1779       s  = TMath::ATan2(   bq*( bq1*(dx*px +dy*py) + a ) , cS )/bq;
1780       ds = TMath::ATan2(TMath::Sqrt(sa),ca)/bq;
1781     } else {
1782       s = ( (dx*px + dy*py) + (py*px1-px*py1)/bq1)/p2;
1783       ds = s*s - (d2-2*(px1*dy-py1*dx)/bq1)/p2; 
1784       if( ds<0 ) ds = 0;
1785       ds = TMath::Sqrt(ds);   
1786     }
1787     
1788     if( TMath::Abs(bq1)>1.e-8){
1789       s1 = TMath::ATan2( -bq1*( bq*(dx*px1+dy*py1) + a), cS1 )/bq1;
1790       ds1 = TMath::ATan2(TMath::Sqrt(sa1),ca1)/bq1;  
1791     } else {
1792       s1 = (-(dx*px1 + dy*py1) + (py*px1-px*py1)/bq)/p21;
1793       ds1 = s1*s1 - (d2+2*(px*dy-py*dx)/bq)/p21; 
1794       if( ds1<0 ) ds1 = 0;
1795       ds1 = TMath::Sqrt(ds1);
1796     }
1797   }
1798
1799   Double_t ss[2], ss1[2], g[2][5],g1[2][5];
1800   
1801   ss[0] = s + ds;
1802   ss[1] = s - ds;
1803   ss1[0] = s1 + ds1;
1804   ss1[1] = s1 - ds1;
1805   for( Int_t i=0; i<2; i++){
1806     Double_t bs = bq*ss[i];
1807     Double_t c = TMath::Cos(bs), sss = TMath::Sin(bs);
1808     Double_t cB,sB;
1809     if( TMath::Abs(bq)>1.e-8){
1810       cB= (1-c)/bq;     
1811       sB= sss/bq;  
1812     }else{
1813       const Double_t kOvSqr6 = 1./TMath::Sqrt(6.);
1814       sB = (1.-bs*kOvSqr6)*(1.+bs*kOvSqr6)*ss[i];
1815       cB = .5*sB*bs;
1816     }
1817     g[i][0] = fP[0] + sB*px + cB*py;
1818     g[i][1] = fP[1] - cB*px + sB*py;
1819     g[i][2] = fP[2] + ss[i]*pz;
1820     g[i][3] =       + c*px + sss*py;
1821     g[i][4] =       - sss*px + c*py;
1822
1823     bs = bq1*ss1[i];  
1824     c =  TMath::Cos(bs); sss = TMath::Sin(bs);
1825     if( TMath::Abs(bq1)>1.e-8){
1826       cB= (1-c)/bq1;   
1827       sB= sss/bq1;  
1828     }else{
1829       const Double_t kOvSqr6 = 1./TMath::Sqrt(6.);
1830       sB = (1.-bs*kOvSqr6)*(1.+bs*kOvSqr6)*ss1[i];
1831       cB = .5*sB*bs;
1832     }
1833       
1834     g1[i][0] = p.fP[0] + sB*px1 + cB*py1;
1835     g1[i][1] = p.fP[1] - cB*px1 + sB*py1;
1836     g1[i][2] = p.fP[2] + ss[i]*pz1;
1837     g1[i][3] =         + c*px1 + sss*py1;
1838     g1[i][4] =         - sss*px1 + c*py1;
1839   }
1840
1841   Int_t i=0, i1=0;
1842   
1843   Double_t dMin = 1.e10;
1844   for( Int_t j=0; j<2; j++){
1845     for( Int_t j1=0; j1<2; j1++){
1846       Double_t xx = g[j][0]-g1[j1][0];
1847       Double_t yy = g[j][1]-g1[j1][1];
1848       Double_t zz = g[j][2]-g1[j1][2];
1849       Double_t d = xx*xx + yy*yy + zz*zz;
1850       if( d<dMin ){
1851         dMin = d;
1852         i = j;
1853         i1 = j1;
1854       }
1855     }
1856   }  
1857
1858   DS = ss[i];
1859   DS1 = ss1[i1];
1860   if(0){
1861     Double_t x= g[i][0], y= g[i][1], z= g[i][2], ppx= g[i][3], ppy= g[i][4];  
1862     Double_t x1=g1[i1][0], y1= g1[i1][1], z1= g1[i1][2], ppx1= g1[i1][3], ppy1= g1[i1][4];  
1863     Double_t dx = x1-x;
1864     Double_t dy = y1-y;
1865     Double_t dz = z1-z;
1866     Double_t a = ppx*ppx1 + ppy*ppy1 + pz*pz1;
1867     Double_t b = dx*ppx1 + dy*ppy1 + dz*pz1;
1868     Double_t c = dx*ppx  + dy*ppy  + dz*pz ;
1869     Double_t pp2 = ppx*ppx + ppy*ppy + pz*pz;
1870     Double_t pp21= ppx1*ppx1 + ppy1*ppy1 + pz1*pz1;
1871     Double_t det = pp2*pp21 - a*a;    
1872     if( TMath::Abs(det)>1.e-8 ){
1873       DS+=(a*b-pp21*c)/det;
1874       DS1+=(a*c-pp2*b)/det;
1875     }
1876   }
1877 }
1878
1879
1880
1881 void AliKFParticleBase::TransportCBM( Double_t dS, 
1882                                  Double_t P[], Double_t C[] ) const
1883 {  
1884   //* Transport the particle on dS, output to P[],C[], for CBM field
1885  
1886   if( fQ==0 ){
1887     TransportLine( dS, P, C );
1888     return;
1889   }
1890
1891   const Double_t kCLight = 0.000299792458;
1892
1893   Double_t c = fQ*kCLight;
1894
1895   // construct coefficients 
1896
1897   Double_t 
1898     px   = fP[3],
1899     py   = fP[4],
1900     pz   = fP[5];
1901       
1902   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;
1903
1904   { // get field integrals
1905
1906     Double_t fld[3][3];   
1907     Double_t p0[3], p1[3], p2[3];
1908
1909     // line track approximation
1910
1911     p0[0] = fP[0];
1912     p0[1] = fP[1];
1913     p0[2] = fP[2];
1914   
1915     p2[0] = fP[0] + px*dS;
1916     p2[1] = fP[1] + py*dS;
1917     p2[2] = fP[2] + pz*dS;
1918   
1919     p1[0] = 0.5*(p0[0]+p2[0]);
1920     p1[1] = 0.5*(p0[1]+p2[1]);
1921     p1[2] = 0.5*(p0[2]+p2[2]);
1922
1923     // first order track approximation
1924     {
1925       GetFieldValue( p0, fld[0] );
1926       GetFieldValue( p1, fld[1] );
1927       GetFieldValue( p2, fld[2] );
1928
1929       Double_t ssy1 = ( 7*fld[0][1] + 6*fld[1][1]-fld[2][1] )*c*dS*dS/96.;
1930       Double_t ssy2 = (   fld[0][1] + 2*fld[1][1]         )*c*dS*dS/6.;
1931
1932       p1[0] -= ssy1*pz;
1933       p1[2] += ssy1*px;
1934       p2[0] -= ssy2*pz;
1935       p2[2] += ssy2*px;   
1936     }
1937
1938     GetFieldValue( p0, fld[0] );
1939     GetFieldValue( p1, fld[1] );
1940     GetFieldValue( p2, fld[2] );
1941     
1942     sx = c*( fld[0][0] + 4*fld[1][0] + fld[2][0] )*dS/6.;
1943     sy = c*( fld[0][1] + 4*fld[1][1] + fld[2][1] )*dS/6.;
1944     sz = c*( fld[0][2] + 4*fld[1][2] + fld[2][2] )*dS/6.;
1945
1946     ssx = c*( fld[0][0] + 2*fld[1][0])*dS*dS/6.;
1947     ssy = c*( fld[0][1] + 2*fld[1][1])*dS*dS/6.;
1948     ssz = c*( fld[0][2] + 2*fld[1][2])*dS*dS/6.;
1949
1950     Double_t c2[3][3]    =   { {  5, -4, -1},{  44,  80,  -4},{ 11, 44, 5} }; // /=360.    
1951     Double_t cc2[3][3]    =   { { 38,  8, -4},{ 148, 208, -20},{  3, 36, 3} }; // /=2520.
1952     for(Int_t n=0; n<3; n++)
1953       for(Int_t m=0; m<3; m++) 
1954         {
1955           syz += c2[n][m]*fld[n][1]*fld[m][2];
1956           ssyz += cc2[n][m]*fld[n][1]*fld[m][2];
1957         }
1958  
1959     syz  *= c*c*dS*dS/360.;
1960     ssyz  *= c*c*dS*dS*dS/2520.;
1961     
1962     syy  = c*( fld[0][1] + 4*fld[1][1] + fld[2][1] )*dS;
1963     syyy = syy*syy*syy / 1296;
1964     syy  = syy*syy/72;
1965
1966     ssyy = ( fld[0][1]*( 38*fld[0][1] + 156*fld[1][1]  -   fld[2][1] )+
1967             fld[1][1]*(              208*fld[1][1]  +16*fld[2][1] )+
1968             fld[2][1]*(                             3*fld[2][1] )  
1969             )*dS*dS*dS*c*c/2520.;
1970     ssyyy = 
1971       (
1972        fld[0][1]*( fld[0][1]*( 85*fld[0][1] + 526*fld[1][1]  - 7*fld[2][1] )+
1973                  fld[1][1]*(             1376*fld[1][1]  +84*fld[2][1] )+
1974                  fld[2][1]*(                            19*fld[2][1] )  )+
1975        fld[1][1]*( fld[1][1]*(             1376*fld[1][1] +256*fld[2][1] )+
1976                  fld[2][1]*(                            62*fld[2][1] )  )+
1977        fld[2][1]*fld[2][1]  *(                             3*fld[2][1] )       
1978        )*dS*dS*dS*dS*c*c*c/90720.;    
1979  
1980   }
1981
1982   Double_t mJ[8][8];
1983   for( Int_t i=0; i<8; i++ ) for( Int_t j=0; j<8; j++) mJ[i][j]=0;
1984
1985   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;
1986   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;
1987   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;
1988   
1989   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;
1990   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;
1991   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;
1992   mJ[6][6] = mJ[7][7] = 1;
1993   
1994   P[0] = fP[0] + mJ[0][3]*px + mJ[0][4]*py + mJ[0][5]*pz;
1995   P[1] = fP[1] + mJ[1][3]*px + mJ[1][4]*py + mJ[1][5]*pz;
1996   P[2] = fP[2] + mJ[2][3]*px + mJ[2][4]*py + mJ[2][5]*pz;
1997   P[3] =        mJ[3][3]*px + mJ[3][4]*py + mJ[3][5]*pz;
1998   P[4] =        mJ[4][3]*px + mJ[4][4]*py + mJ[4][5]*pz;
1999   P[5] =        mJ[5][3]*px + mJ[5][4]*py + mJ[5][5]*pz;
2000   P[6] = fP[6];
2001   P[7] = fP[7];
2002
2003   MultQSQt( mJ[0], fC, C);
2004
2005 }
2006
2007
2008 void AliKFParticleBase::TransportBz( Double_t b, Double_t t,
2009                                      Double_t p[], Double_t e[] ) const 
2010
2011   //* Transport the particle on dS, output to P[],C[], for Bz field
2012  
2013   const Double_t kCLight = 0.000299792458;
2014   b = b*fQ*kCLight;
2015   Double_t bs= b*t;
2016   Double_t s = TMath::Sin(bs), c = TMath::Cos(bs);
2017   Double_t sB, cB;
2018   if( TMath::Abs(bs)>1.e-10){
2019     sB= s/b;
2020     cB= (1-c)/b;
2021   }else{
2022     const Double_t kOvSqr6 = 1./TMath::Sqrt(6.);
2023     sB = (1.-bs*kOvSqr6)*(1.+bs*kOvSqr6)*t;
2024     cB = .5*sB*bs;
2025   }
2026   
2027   Double_t px = fP[3];
2028   Double_t py = fP[4];
2029   Double_t pz = fP[5];
2030   
2031   p[0] = fP[0] + sB*px + cB*py;
2032   p[1] = fP[1] - cB*px + sB*py;
2033   p[2] = fP[2] +  t*pz;
2034   p[3] =          c*px + s*py;
2035   p[4] =         -s*px + c*py;
2036   p[5] = fP[5];
2037   p[6] = fP[6];
2038   p[7] = fP[7];
2039
2040   /* 
2041   Double_t mJ[8][8] = { {1,0,0,   sB, cB,  0, 0, 0 },
2042                         {0,1,0,  -cB, sB,  0, 0, 0 },
2043                         {0,0,1,    0,  0,  t, 0, 0 },
2044                         {0,0,0,    c,  s,  0, 0, 0 },
2045                         {0,0,0,   -s,  c,  0, 0, 0 },
2046                         {0,0,0,    0,  0,  1, 0, 0 },
2047                         {0,0,0,    0,  0,  0, 1, 0 },
2048                         {0,0,0,    0,  0,  0, 0, 1 }  };
2049   Double_t mA[8][8];
2050   for( Int_t k=0,i=0; i<8; i++)
2051     for( Int_t j=0; j<=i; j++, k++ ) mA[i][j] = mA[j][i] = fC[k]; 
2052
2053   Double_t mJC[8][8];
2054   for( Int_t i=0; i<8; i++ )
2055     for( Int_t j=0; j<8; j++ ){
2056       mJC[i][j]=0;
2057       for( Int_t k=0; k<8; k++ ) mJC[i][j]+=mJ[i][k]*mA[k][j];
2058     }
2059   
2060   for( Int_t k=0,i=0; i<8; i++)
2061     for( Int_t j=0; j<=i; j++, k++ ){
2062       e[k] = 0;
2063       for( Int_t l=0; l<8; l++ ) e[k]+=mJC[i][l]*mJ[j][l];
2064     }
2065   
2066   return;
2067   */
2068
2069   Double_t 
2070     c6=fC[6], c7=fC[7], c8=fC[8], c17=fC[17], c18=fC[18],
2071     c24 = fC[24], c31 = fC[31];
2072
2073   Double_t 
2074     cBC13 = cB*fC[13],
2075     mJC13 = c7 - cB*fC[9] + sB*fC[13],
2076     mJC14 = fC[11] - cBC13 + sB*fC[14],
2077     mJC23 = c8 + t*c18,
2078     mJC24 = fC[12] + t*fC[19],
2079     mJC33 = c*fC[9] + s*fC[13],
2080     mJC34 = c*fC[13] + s*fC[14],
2081     mJC43 = -s*fC[9] + c*fC[13],
2082     mJC44 = -s*fC[13] + c*fC[14];
2083
2084
2085   e[0]= fC[0] + 2*(sB*c6 + cB*fC[10]) + (sB*fC[9] + 2*cBC13)*sB + cB*cB*fC[14];
2086   e[1]= fC[1] - cB*c6 + sB*fC[10] + mJC13*sB + mJC14*cB;
2087   e[2]= fC[2] - cB*c7 + sB*fC[11] - mJC13*cB + mJC14*sB;
2088   e[3]= fC[3] + t*fC[15] + mJC23*sB + mJC24*cB;
2089   e[4]= fC[4] + t*fC[16] - mJC23*cB + mJC24*sB;
2090
2091   e[15]= fC[15] + c18*sB + fC[19]*cB;
2092   e[16]= fC[16] - c18*cB + fC[19]*sB;
2093   e[17]= c17 + fC[20]*t;
2094   e[18]= c18*c + fC[19]*s;
2095   e[19]= -c18*s + fC[19]*c;
2096
2097   e[5]= fC[5] + (c17 + e[17] )*t;
2098
2099   e[6]= c*c6 + s*fC[10] + mJC33*sB + mJC34*cB;
2100   e[7]= c*c7 + s*fC[11] - mJC33*cB + mJC34*sB;
2101   e[8]= c*c8 + s*fC[12] + e[18]*t;
2102   e[9]= mJC33*c + mJC34*s;
2103   e[10]= -s*c6 + c*fC[10] + mJC43*sB + mJC44*cB;
2104
2105     
2106   e[11]= -s*c7 + c*fC[11] - mJC43*cB + mJC44*sB;
2107   e[12]= -s*c8 + c*fC[12] + e[19]*t;
2108   e[13]= mJC43*c + mJC44*s;
2109   e[14]= -mJC43*s + mJC44*c;
2110   e[20]= fC[20];
2111   e[21]= fC[21] + fC[25]*cB + c24*sB;
2112   e[22]= fC[22] - c24*cB + fC[25]*sB;
2113   e[23]= fC[23] + fC[26]*t;
2114   e[24]= c*c24 + s*fC[25];
2115   e[25]= c*fC[25] - c24*s;
2116   e[26]= fC[26];
2117   e[27]= fC[27];
2118   e[28]= fC[28] + fC[32]*cB + c31*sB;
2119   e[29]= fC[29] - c31*cB + fC[32]*sB;
2120   e[30]= fC[30] + fC[33]*t;
2121   e[31]= c*c31 + s*fC[32];
2122   e[32]= c*fC[32] - s*c31;
2123   e[33]= fC[33];
2124   e[34]= fC[34];
2125   e[35]= fC[35];     
2126 }
2127
2128
2129 Double_t AliKFParticleBase::GetDistanceFromVertex( const AliKFParticleBase &Vtx ) const
2130 {
2131   //* Calculate distance from vertex [cm]
2132
2133   return GetDistanceFromVertex( Vtx.fP );
2134 }
2135
2136 Double_t AliKFParticleBase::GetDistanceFromVertex( const Double_t vtx[] ) const
2137 {
2138   //* Calculate distance from vertex [cm]
2139
2140   Double_t mP[8], mC[36];  
2141   Transport( GetDStoPoint(vtx), mP, mC );
2142   Double_t d[3]={ vtx[0]-mP[0], vtx[1]-mP[1], vtx[2]-mP[2]};
2143   return TMath::Sqrt( d[0]*d[0]+d[1]*d[1]+d[2]*d[2] );
2144 }
2145
2146 Double_t AliKFParticleBase::GetDistanceFromParticle( const AliKFParticleBase &p ) 
2147   const
2148
2149   //* Calculate distance to other particle [cm]
2150
2151   Double_t dS, dS1;
2152   GetDStoParticle( p, dS, dS1 );   
2153   Double_t mP[8], mC[36], mP1[8], mC1[36];
2154   Transport( dS, mP, mC ); 
2155   p.Transport( dS1, mP1, mC1 ); 
2156   Double_t dx = mP[0]-mP1[0]; 
2157   Double_t dy = mP[1]-mP1[1]; 
2158   Double_t dz = mP[2]-mP1[2]; 
2159   dz = 0;
2160   return TMath::Sqrt(dx*dx+dy*dy+dz*dz);
2161 }
2162
2163 Double_t AliKFParticleBase::GetDeviationFromVertex( const AliKFParticleBase &Vtx ) const
2164 {
2165   //* Calculate sqrt(Chi2/ndf) deviation from vertex
2166
2167   return GetDeviationFromVertex( Vtx.fP, Vtx.fC );
2168 }
2169
2170
2171 Double_t AliKFParticleBase::GetDeviationFromVertex( const Double_t v[], const Double_t Cv[] ) const
2172 {
2173   //* Calculate sqrt(Chi2/ndf) deviation from vertex
2174   //* v = [xyz], Cv=[Cxx,Cxy,Cyy,Cxz,Cyz,Czz]-covariance matrix
2175
2176   Double_t mP[8];
2177   Double_t mC[36];
2178   
2179   Transport( GetDStoPoint(v), mP, mC );  
2180
2181   Double_t d[3]={ v[0]-mP[0], v[1]-mP[1], v[2]-mP[2]};
2182
2183   Double_t sigmaS = .1+10.*TMath::Sqrt( (d[0]*d[0]+d[1]*d[1]+d[2]*d[2])/
2184                                  (mP[3]*mP[3]+mP[4]*mP[4]+mP[5]*mP[5])  );
2185
2186    
2187   Double_t h[3] = { mP[3]*sigmaS, mP[4]*sigmaS, mP[5]*sigmaS };       
2188   
2189   Double_t mSi[6] = 
2190     { mC[0] +h[0]*h[0], 
2191       mC[1] +h[1]*h[0], mC[2] +h[1]*h[1], 
2192       mC[3] +h[2]*h[0], mC[4] +h[2]*h[1], mC[5] +h[2]*h[2] };
2193
2194   if( Cv ){
2195     mSi[0]+=Cv[0];
2196     mSi[1]+=Cv[1];
2197     mSi[2]+=Cv[2];
2198     mSi[3]+=Cv[3];
2199     mSi[4]+=Cv[4];
2200     mSi[5]+=Cv[5];
2201   }
2202   
2203   Double_t mS[6];
2204
2205   mS[0] = mSi[2]*mSi[5] - mSi[4]*mSi[4];
2206   mS[1] = mSi[3]*mSi[4] - mSi[1]*mSi[5];
2207   mS[2] = mSi[0]*mSi[5] - mSi[3]*mSi[3];
2208   mS[3] = mSi[1]*mSi[4] - mSi[2]*mSi[3];
2209   mS[4] = mSi[1]*mSi[3] - mSi[0]*mSi[4];
2210   mS[5] = mSi[0]*mSi[2] - mSi[1]*mSi[1];         
2211       
2212   Double_t s = ( mSi[0]*mS[0] + mSi[1]*mS[1] + mSi[3]*mS[3] );
2213   s = ( s > 1.E-20 )  ?1./s :0;   
2214
2215   return TMath::Sqrt( TMath::Abs(s*( ( mS[0]*d[0] + mS[1]*d[1] + mS[3]*d[2])*d[0]
2216                    +(mS[1]*d[0] + mS[2]*d[1] + mS[4]*d[2])*d[1]
2217                    +(mS[3]*d[0] + mS[4]*d[1] + mS[5]*d[2])*d[2] ))/2);
2218 }
2219
2220
2221 Double_t AliKFParticleBase::GetDeviationFromParticle( const AliKFParticleBase &p ) 
2222   const
2223
2224   //* Calculate sqrt(Chi2/ndf) deviation from other particle
2225
2226   Double_t dS, dS1;
2227   GetDStoParticle( p, dS, dS1 );   
2228   Double_t mP1[8], mC1[36];
2229   p.Transport( dS1, mP1, mC1 ); 
2230
2231   Double_t d[3]={ fP[0]-mP1[0], fP[1]-mP1[1], fP[2]-mP1[2]};
2232
2233   Double_t sigmaS = .1+10.*TMath::Sqrt( (d[0]*d[0]+d[1]*d[1]+d[2]*d[2])/
2234                                         (mP1[3]*mP1[3]+mP1[4]*mP1[4]+mP1[5]*mP1[5])  );
2235
2236   Double_t h[3] = { mP1[3]*sigmaS, mP1[4]*sigmaS, mP1[5]*sigmaS };       
2237   
2238   mC1[0] +=h[0]*h[0];
2239   mC1[1] +=h[1]*h[0]; 
2240   mC1[2] +=h[1]*h[1]; 
2241   mC1[3] +=h[2]*h[0]; 
2242   mC1[4] +=h[2]*h[1];
2243   mC1[5] +=h[2]*h[2];
2244
2245   return GetDeviationFromVertex( mP1, mC1 )*TMath::Sqrt(2./1.);
2246 }
2247
2248
2249
2250 void AliKFParticleBase::SubtractFromVertex(  AliKFParticleBase &Vtx ) const
2251 {    
2252   //* Subtract the particle from the vertex  
2253
2254   Double_t fld[3];  
2255   {
2256     GetFieldValue( Vtx.fP, fld );
2257     const Double_t kCLight =  0.000299792458;
2258     fld[0]*=kCLight; fld[1]*=kCLight; fld[2]*=kCLight;
2259   }
2260
2261   Double_t m[8];
2262   Double_t mCm[36];
2263
2264   if( Vtx.fIsLinearized ){
2265     GetMeasurement( Vtx.fVtxGuess, m, mCm );
2266   } else {
2267     GetMeasurement( Vtx.fP, m, mCm );
2268   }
2269   
2270   Double_t mV[6];
2271     
2272   mV[ 0] = mCm[ 0];
2273   mV[ 1] = mCm[ 1];
2274   mV[ 2] = mCm[ 2];
2275   mV[ 3] = mCm[ 3];
2276   mV[ 4] = mCm[ 4];
2277   mV[ 5] = mCm[ 5];
2278      
2279   //* 
2280             
2281   Double_t mS[6];
2282   {
2283     Double_t mSi[6] = { mV[0]-Vtx.fC[0], 
2284                         mV[1]-Vtx.fC[1], mV[2]-Vtx.fC[2], 
2285                         mV[3]-Vtx.fC[3], mV[4]-Vtx.fC[4], mV[5]-Vtx.fC[5] };
2286     
2287     mS[0] = mSi[2]*mSi[5] - mSi[4]*mSi[4];
2288     mS[1] = mSi[3]*mSi[4] - mSi[1]*mSi[5];
2289     mS[2] = mSi[0]*mSi[5] - mSi[3]*mSi[3];
2290     mS[3] = mSi[1]*mSi[4] - mSi[2]*mSi[3];
2291     mS[4] = mSi[1]*mSi[3] - mSi[0]*mSi[4];
2292     mS[5] = mSi[0]*mSi[2] - mSi[1]*mSi[1];       
2293     
2294     Double_t s = ( mSi[0]*mS[0] + mSi[1]*mS[1] + mSi[3]*mS[3] );
2295     s = ( s > 1.E-20 )  ?1./s :0;         
2296     mS[0]*=s;
2297     mS[1]*=s;
2298     mS[2]*=s;
2299     mS[3]*=s;
2300     mS[4]*=s;
2301     mS[5]*=s;
2302   }
2303     
2304   //* Residual (measured - estimated)
2305     
2306   Double_t zeta[3] = { m[0]-Vtx.fP[0], m[1]-Vtx.fP[1], m[2]-Vtx.fP[2] };
2307         
2308   //* mCHt = mCH' - D'
2309     
2310   Double_t mCHt0[3], mCHt1[3], mCHt2[3];
2311     
2312   mCHt0[0]=Vtx.fC[ 0] ;      mCHt1[0]=Vtx.fC[ 1] ;      mCHt2[0]=Vtx.fC[ 3] ;
2313   mCHt0[1]=Vtx.fC[ 1] ;      mCHt1[1]=Vtx.fC[ 2] ;      mCHt2[1]=Vtx.fC[ 4] ;
2314   mCHt0[2]=Vtx.fC[ 3] ;      mCHt1[2]=Vtx.fC[ 4] ;      mCHt2[2]=Vtx.fC[ 5] ;
2315   
2316   //* Kalman gain K = mCH'*S
2317     
2318   Double_t k0[3], k1[3], k2[3];
2319     
2320   for(Int_t i=0;i<3;++i){
2321     k0[i] = mCHt0[i]*mS[0] + mCHt1[i]*mS[1] + mCHt2[i]*mS[3];
2322     k1[i] = mCHt0[i]*mS[1] + mCHt1[i]*mS[2] + mCHt2[i]*mS[4];
2323     k2[i] = mCHt0[i]*mS[3] + mCHt1[i]*mS[4] + mCHt2[i]*mS[5];
2324   }
2325     
2326   //* New estimation of the vertex position r += K*zeta
2327     
2328   Double_t dChi2 = -(mS[0]*zeta[0] + mS[1]*zeta[1] + mS[3]*zeta[2])*zeta[0]
2329     +      (mS[1]*zeta[0] + mS[2]*zeta[1] + mS[4]*zeta[2])*zeta[1]
2330     +      (mS[3]*zeta[0] + mS[4]*zeta[1] + mS[5]*zeta[2])*zeta[2];
2331
2332   if( Vtx.fChi2 - dChi2 < 0 ) return;
2333
2334   for(Int_t i=0;i<3;++i) 
2335     Vtx.fP[i] -= k0[i]*zeta[0] + k1[i]*zeta[1] + k2[i]*zeta[2];       
2336     
2337   //* New covariance matrix C -= K*(mCH')'
2338     
2339   for(Int_t i=0, k=0;i<3;++i){
2340     for(Int_t j=0;j<=i;++j,++k) 
2341       Vtx.fC[k] += k0[i]*mCHt0[j] + k1[i]*mCHt1[j] + k2[i]*mCHt2[j];
2342   }
2343     
2344   //* Calculate Chi^2 
2345
2346   Vtx.fNDF  -= 2;
2347   Vtx.fChi2 -= dChi2;
2348 }
2349
2350
2351
2352 void AliKFParticleBase::TransportLine( Double_t dS, 
2353                                        Double_t P[], Double_t C[] ) const 
2354 {
2355   //* Transport the particle as a straight line
2356
2357   P[0] = fP[0] + dS*fP[3];
2358   P[1] = fP[1] + dS*fP[4];
2359   P[2] = fP[2] + dS*fP[5];
2360   P[3] = fP[3];
2361   P[4] = fP[4];
2362   P[5] = fP[5];
2363   P[6] = fP[6];
2364   P[7] = fP[7];
2365  
2366   Double_t c6  = fC[ 6] + dS*fC[ 9];
2367   Double_t c11 = fC[11] + dS*fC[14];
2368   Double_t c17 = fC[17] + dS*fC[20];
2369   Double_t sc13 = dS*fC[13];
2370   Double_t sc18 = dS*fC[18];
2371   Double_t sc19 = dS*fC[19];
2372
2373   C[ 0] = fC[ 0] + dS*( fC[ 6] + c6  );
2374   C[ 2] = fC[ 2] + dS*( fC[11] + c11 );
2375   C[ 5] = fC[ 5] + dS*( fC[17] + c17 );
2376
2377   C[ 7] = fC[ 7] + sc13;
2378   C[ 8] = fC[ 8] + sc18;
2379   C[ 9] = fC[ 9];
2380
2381   C[12] = fC[12] + sc19;
2382
2383   C[ 1] = fC[ 1] + dS*( fC[10] + C[ 7] );
2384   C[ 3] = fC[ 3] + dS*( fC[15] + C[ 8] );
2385   C[ 4] = fC[ 4] + dS*( fC[16] + C[12] ); 
2386   C[ 6] = c6;
2387
2388   C[10] = fC[10] + sc13;
2389   C[11] = c11;
2390
2391   C[13] = fC[13];
2392   C[14] = fC[14];
2393   C[15] = fC[15] + sc18;
2394   C[16] = fC[16] + sc19;
2395   C[17] = c17;
2396   
2397   C[18] = fC[18];
2398   C[19] = fC[19];
2399   C[20] = fC[20];
2400   C[21] = fC[21] + dS*fC[24];
2401   C[22] = fC[22] + dS*fC[25];
2402   C[23] = fC[23] + dS*fC[26];
2403
2404   C[24] = fC[24];
2405   C[25] = fC[25];
2406   C[26] = fC[26];
2407   C[27] = fC[27];
2408   C[28] = fC[28] + dS*fC[31];
2409   C[29] = fC[29] + dS*fC[32];
2410   C[30] = fC[30] + dS*fC[33];
2411
2412   C[31] = fC[31];
2413   C[32] = fC[32];
2414   C[33] = fC[33];
2415   C[34] = fC[34];
2416   C[35] = fC[35]; 
2417 }
2418
2419
2420 void AliKFParticleBase::ConstructGammaBz( const AliKFParticleBase &daughter1,
2421                                           const AliKFParticleBase &daughter2, double Bz  )
2422
2423   //* Create gamma
2424   
2425   const AliKFParticleBase *daughters[2] = { &daughter1, &daughter2};
2426
2427   double v0[3];
2428   
2429   if( !fIsLinearized ){
2430     Double_t ds, ds1;
2431     Double_t m[8];
2432     Double_t mCd[36];       
2433     daughter1.GetDStoParticle(daughter2, ds, ds1);      
2434     daughter1.Transport( ds, m, mCd );
2435     fP[0] = m[0];
2436     fP[1] = m[1];
2437     fP[2] = m[2];
2438     daughter2.Transport( ds1, m, mCd );
2439     fP[0] = .5*( fP[0] + m[0] );
2440     fP[1] = .5*( fP[1] + m[1] );
2441     fP[2] = .5*( fP[2] + m[2] );
2442   } else {
2443     fP[0] = fVtxGuess[0];
2444     fP[1] = fVtxGuess[1];
2445     fP[2] = fVtxGuess[2];
2446   }
2447
2448   double daughterP[2][8], daughterC[2][36];
2449   double vtxMom[2][3];
2450
2451   int nIter = fIsLinearized ?1 :2;
2452
2453   for( int iter=0; iter<nIter; iter++){
2454
2455     v0[0] = fP[0];
2456     v0[1] = fP[1];
2457     v0[2] = fP[2];
2458     
2459     fAtProductionVertex = 0;
2460     fSFromDecay = 0;
2461     fP[0] = v0[0];
2462     fP[1] = v0[1];
2463     fP[2] = v0[2];
2464     fP[3] = 0;
2465     fP[4] = 0;
2466     fP[5] = 0;
2467     fP[6] = 0;
2468     fP[7] = 0;
2469
2470   
2471     // fit daughters to the vertex guess  
2472     
2473     {  
2474       for( int id=0; id<2; id++ ){
2475         
2476         double *p = daughterP[id];
2477         double *mC = daughterC[id];
2478         
2479         daughters[id]->GetMeasurement( v0, p, mC );
2480         
2481         Double_t mAi[6];
2482         InvertSym3(mC, mAi );
2483         
2484         Double_t mB[3][3];
2485         
2486         mB[0][0] = mC[ 6]*mAi[0] + mC[ 7]*mAi[1] + mC[ 8]*mAi[3];
2487         mB[0][1] = mC[ 6]*mAi[1] + mC[ 7]*mAi[2] + mC[ 8]*mAi[4];
2488         mB[0][2] = mC[ 6]*mAi[3] + mC[ 7]*mAi[4] + mC[ 8]*mAi[5];
2489         
2490         mB[1][0] = mC[10]*mAi[0] + mC[11]*mAi[1] + mC[12]*mAi[3];
2491         mB[1][1] = mC[10]*mAi[1] + mC[11]*mAi[2] + mC[12]*mAi[4];
2492         mB[1][2] = mC[10]*mAi[3] + mC[11]*mAi[4] + mC[12]*mAi[5];
2493         
2494         mB[2][0] = mC[15]*mAi[0] + mC[16]*mAi[1] + mC[17]*mAi[3];
2495         mB[2][1] = mC[15]*mAi[1] + mC[16]*mAi[2] + mC[17]*mAi[4];
2496         mB[2][2] = mC[15]*mAi[3] + mC[16]*mAi[4] + mC[17]*mAi[5];
2497         
2498         Double_t z[3] = { v0[0]-p[0], v0[1]-p[1], v0[2]-p[2] };
2499         
2500         vtxMom[id][0] = p[3] + mB[0][0]*z[0] + mB[0][1]*z[1] + mB[0][2]*z[2];
2501         vtxMom[id][1] = p[4] + mB[1][0]*z[0] + mB[1][1]*z[1] + mB[1][2]*z[2];
2502         vtxMom[id][2] = p[5] + mB[2][0]*z[0] + mB[2][1]*z[1] + mB[2][2]*z[2];
2503         
2504         daughters[id]->Transport( daughters[id]->GetDStoPoint(v0), p, mC );
2505       
2506       }      
2507       
2508     } // fit daughters to guess
2509
2510     
2511     // fit new vertex
2512     {
2513       
2514       double mpx0 =  vtxMom[0][0]+vtxMom[1][0];
2515       double mpy0 =  vtxMom[0][1]+vtxMom[1][1];
2516       double mpt0 = TMath::Sqrt(mpx0*mpx0 + mpy0*mpy0);
2517       // double a0 = TMath::ATan2(mpy0,mpx0);
2518       
2519       double ca0 = mpx0/mpt0;
2520       double sa0 = mpy0/mpt0;
2521       double r[3] = { v0[0], v0[1], v0[2] };
2522       double mC[3][3] = {{10000., 0 ,   0  },
2523                          {0,  10000.,   0  },
2524                          {0,     0, 10000. } };
2525       double chi2=0;
2526       
2527       for( int id=0; id<2; id++ ){              
2528         const Double_t kCLight = 0.000299792458;
2529         Double_t q = Bz*daughters[id]->GetQ()*kCLight;
2530         Double_t px0 = vtxMom[id][0];
2531         Double_t py0 = vtxMom[id][1];
2532         Double_t pz0 = vtxMom[id][2];
2533         Double_t pt0 = TMath::Sqrt(px0*px0+py0*py0);
2534         Double_t mG[3][6], mB[3], mH[3][3];
2535         // r = {vx,vy,vz};
2536         // m = {x,y,z,Px,Py,Pz};
2537         // V = daughter.C
2538         // G*m + B = H*r;
2539         // q*x + Py - q*vx - sin(a)*Pt = 0
2540         // q*y - Px - q*vy + cos(a)*Pt = 0
2541         // (Px*cos(a) + Py*sin(a) ) (vz -z) - Pz( cos(a)*(vx-x) + sin(a)*(vy-y)) = 0
2542         
2543         mG[0][0] = q;
2544         mG[0][1] = 0;
2545         mG[0][2] = 0;
2546         mG[0][3] =   -sa0*px0/pt0;
2547         mG[0][4] = 1 -sa0*py0/pt0;
2548         mG[0][5] = 0;   
2549         mH[0][0] = q;
2550         mH[0][1] = 0;
2551         mH[0][2] = 0;      
2552         mB[0] = py0 - sa0*pt0 - mG[0][3]*px0 - mG[0][4]*py0 ;
2553         
2554         // q*y - Px - q*vy + cos(a)*Pt = 0
2555         
2556         mG[1][0] = 0;
2557         mG[1][1] = q;
2558         mG[1][2] = 0;
2559         mG[1][3] = -1 + ca0*px0/pt0;
2560         mG[1][4] =    + ca0*py0/pt0;
2561         mG[1][5] = 0;      
2562         mH[1][0] = 0;
2563         mH[1][1] = q;
2564         mH[1][2] = 0;      
2565         mB[1] = -px0 + ca0*pt0 - mG[1][3]*px0 - mG[1][4]*py0 ;
2566         
2567         // (Px*cos(a) + Py*sin(a) ) (z -vz) - Pz( cos(a)*(x-vx) + sin(a)*(y-vy)) = 0
2568       
2569         mG[2][0] = -pz0*ca0;
2570         mG[2][1] = -pz0*sa0;
2571         mG[2][2] =  px0*ca0 + py0*sa0;
2572         mG[2][3] = 0;
2573         mG[2][4] = 0;
2574         mG[2][5] = 0;
2575         
2576         mH[2][0] = mG[2][0];
2577         mH[2][1] = mG[2][1];
2578         mH[2][2] = mG[2][2];
2579         
2580         mB[2] = 0;
2581         
2582         // fit the vertex
2583
2584         // V = GVGt
2585
2586         double mGV[3][6];
2587         double mV[6];
2588         double m[3];
2589         for( int i=0; i<3; i++ ){
2590           m[i] = mB[i];
2591           for( int k=0; k<6; k++ ) m[i]+=mG[i][k]*daughterP[id][k];
2592         }
2593         for( int i=0; i<3; i++ ){
2594           for( int j=0; j<6; j++ ){
2595             mGV[i][j] = 0;
2596             for( int k=0; k<6; k++ ) mGV[i][j]+=mG[i][k]*daughterC[id][ IJ(k,j) ];
2597           }
2598         }
2599         for( int i=0, k=0; i<3; i++ ){
2600           for( int j=0; j<=i; j++,k++ ){
2601             mV[k] = 0;
2602             for( int l=0; l<6; l++ ) mV[k]+=mGV[i][l]*mG[j][l];
2603           }
2604         }
2605         
2606       
2607         //* CHt
2608         
2609         Double_t mCHt[3][3];
2610         Double_t mHCHt[6];
2611         Double_t mHr[3];
2612         for( int i=0; i<3; i++ ){         
2613           mHr[i] = 0;
2614           for( int k=0; k<3; k++ ) mHr[i]+= mH[i][k]*r[k];
2615         }
2616       
2617         for( int i=0; i<3; i++ ){
2618           for( int j=0; j<3; j++){
2619             mCHt[i][j] = 0;
2620             for( int k=0; k<3; k++ ) mCHt[i][j]+= mC[i][k]*mH[j][k];
2621           }
2622         }
2623
2624         for( int i=0, k=0; i<3; i++ ){
2625           for( int j=0; j<=i; j++, k++ ){
2626             mHCHt[k] = 0;
2627             for( int l=0; l<3; l++ ) mHCHt[k]+= mH[i][l]*mCHt[l][j];
2628           }
2629         }
2630       
2631         Double_t mS[6] = { mHCHt[0]+mV[0], 
2632                            mHCHt[1]+mV[1], mHCHt[2]+mV[2], 
2633                            mHCHt[3]+mV[3], mHCHt[4]+mV[4], mHCHt[5]+mV[5]    }; 
2634       
2635
2636         InvertSym3(mS,mS);
2637         
2638         //* Residual (measured - estimated)
2639     
2640         Double_t zeta[3] = { m[0]-mHr[0], m[1]-mHr[1], m[2]-mHr[2] };
2641             
2642         //* Kalman gain K = mCH'*S
2643     
2644         Double_t k[3][3];
2645       
2646         for(Int_t i=0;i<3;++i){
2647           k[i][0] = mCHt[i][0]*mS[0] + mCHt[i][1]*mS[1] + mCHt[i][2]*mS[3];
2648           k[i][1] = mCHt[i][0]*mS[1] + mCHt[i][1]*mS[2] + mCHt[i][2]*mS[4];
2649           k[i][2] = mCHt[i][0]*mS[3] + mCHt[i][1]*mS[4] + mCHt[i][2]*mS[5];
2650         }
2651
2652         //* New estimation of the vertex position r += K*zeta
2653     
2654         for(Int_t i=0;i<3;++i) 
2655           r[i] = r[i] + k[i][0]*zeta[0] + k[i][1]*zeta[1] + k[i][2]*zeta[2];
2656       
2657         //* New covariance matrix C -= K*(mCH')'
2658
2659         for(Int_t i=0;i<3;++i){
2660           for(Int_t j=0;j<=i;++j){
2661             mC[i][j] = mC[i][j] - (k[i][0]*mCHt[j][0] + k[i][1]*mCHt[j][1] + k[i][2]*mCHt[j][2]);
2662             mC[j][i] = mC[i][j];
2663           }
2664         }
2665
2666         
2667         //* Calculate Chi^2 
2668         
2669         chi2 += ( ( mS[0]*zeta[0] + mS[1]*zeta[1] + mS[3]*zeta[2] )*zeta[0]
2670                   +(mS[1]*zeta[0] + mS[2]*zeta[1] + mS[4]*zeta[2] )*zeta[1]
2671                   +(mS[3]*zeta[0] + mS[4]*zeta[1] + mS[5]*zeta[2] )*zeta[2]  );  
2672       }
2673     
2674       // store vertex
2675     
2676       fNDF  = 2;
2677       fChi2 = chi2;
2678       for( int i=0; i<3; i++ ) fP[i] = r[i];
2679       for( int i=0,k=0; i<3; i++ ){
2680         for( int j=0; j<=i; j++,k++ ){
2681           fC[k] = mC[i][j];
2682         }
2683       }
2684     }
2685
2686   } // iterations
2687
2688   // now fit daughters to the vertex
2689   
2690   fQ     =  0;
2691   fSFromDecay = 0;    
2692
2693   for(Int_t i=3;i<8;++i) fP[i]=0.;
2694   for(Int_t i=6;i<35;++i) fC[i]=0.;
2695   fC[35] = 100.;
2696
2697   for( int id=0; id<2; id++ ){
2698
2699     double *p = daughterP[id];
2700     double *mC = daughterC[id];      
2701     daughters[id]->GetMeasurement( v0, p, mC );
2702
2703     const Double_t *m = fP, *mV = fC;
2704     
2705     Double_t mAi[6];
2706     InvertSym3(mC, mAi );
2707
2708     Double_t mB[4][3];
2709
2710     mB[0][0] = mC[ 6]*mAi[0] + mC[ 7]*mAi[1] + mC[ 8]*mAi[3];
2711     mB[0][1] = mC[ 6]*mAi[1] + mC[ 7]*mAi[2] + mC[ 8]*mAi[4];
2712     mB[0][2] = mC[ 6]*mAi[3] + mC[ 7]*mAi[4] + mC[ 8]*mAi[5];
2713     
2714     mB[1][0] = mC[10]*mAi[0] + mC[11]*mAi[1] + mC[12]*mAi[3];
2715     mB[1][1] = mC[10]*mAi[1] + mC[11]*mAi[2] + mC[12]*mAi[4];
2716     mB[1][2] = mC[10]*mAi[3] + mC[11]*mAi[4] + mC[12]*mAi[5];
2717     
2718     mB[2][0] = mC[15]*mAi[0] + mC[16]*mAi[1] + mC[17]*mAi[3];
2719     mB[2][1] = mC[15]*mAi[1] + mC[16]*mAi[2] + mC[17]*mAi[4];
2720     mB[2][2] = mC[15]*mAi[3] + mC[16]*mAi[4] + mC[17]*mAi[5];
2721     
2722     mB[3][0] = mC[21]*mAi[0] + mC[22]*mAi[1] + mC[23]*mAi[3];
2723     mB[3][1] = mC[21]*mAi[1] + mC[22]*mAi[2] + mC[23]*mAi[4];
2724     mB[3][2] = mC[21]*mAi[3] + mC[22]*mAi[4] + mC[23]*mAi[5];    
2725
2726
2727     Double_t z[3] = { m[0]-p[0], m[1]-p[1], m[2]-p[2] };
2728
2729     {
2730       Double_t mAV[6] = { mC[0]-mV[0], mC[1]-mV[1], mC[2]-mV[2], 
2731                           mC[3]-mV[3], mC[4]-mV[4], mC[5]-mV[5] };
2732       
2733       Double_t mAVi[6];
2734       if( !InvertSym3(mAV, mAVi) ){
2735         Double_t dChi2 = ( +(mAVi[0]*z[0] + mAVi[1]*z[1] + mAVi[3]*z[2])*z[0]
2736                            +(mAVi[1]*z[0] + mAVi[2]*z[1] + mAVi[4]*z[2])*z[1]
2737                            +(mAVi[3]*z[0] + mAVi[4]*z[1] + mAVi[5]*z[2])*z[2] );
2738         fChi2+= TMath::Abs( dChi2 );
2739       }
2740       fNDF  += 2;
2741     }
2742
2743     //* Add the daughter momentum to the particle momentum
2744  
2745     fP[3]+= p[3] + mB[0][0]*z[0] + mB[0][1]*z[1] + mB[0][2]*z[2];
2746     fP[4]+= p[4] + mB[1][0]*z[0] + mB[1][1]*z[1] + mB[1][2]*z[2];
2747     fP[5]+= p[5] + mB[2][0]*z[0] + mB[2][1]*z[1] + mB[2][2]*z[2];
2748     fP[6]+= p[6] + mB[3][0]*z[0] + mB[3][1]*z[1] + mB[3][2]*z[2];
2749   
2750     Double_t d0, d1, d2;
2751    
2752     d0= mB[0][0]*mV[0] + mB[0][1]*mV[1] + mB[0][2]*mV[3] - mC[ 6];
2753     d1= mB[0][0]*mV[1] + mB[0][1]*mV[2] + mB[0][2]*mV[4] - mC[ 7];
2754     d2= mB[0][0]*mV[3] + mB[0][1]*mV[4] + mB[0][2]*mV[5] - mC[ 8];
2755
2756     //fC[6]+= mC[ 6] + d0;
2757     //fC[7]+= mC[ 7] + d1;
2758     //fC[8]+= mC[ 8] + d2;
2759     fC[9]+= mC[ 9] + d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
2760
2761     d0= mB[1][0]*mV[0] + mB[1][1]*mV[1] + mB[1][2]*mV[3] - mC[10];
2762     d1= mB[1][0]*mV[1] + mB[1][1]*mV[2] + mB[1][2]*mV[4] - mC[11];
2763     d2= mB[1][0]*mV[3] + mB[1][1]*mV[4] + mB[1][2]*mV[5] - mC[12];
2764
2765     //fC[10]+= mC[10]+ d0;
2766     //fC[11]+= mC[11]+ d1;
2767     //fC[12]+= mC[12]+ d2;
2768     fC[13]+= mC[13]+ d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
2769     fC[14]+= mC[14]+ d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
2770
2771     d0= mB[2][0]*mV[0] + mB[2][1]*mV[1] + mB[2][2]*mV[3] - mC[15];
2772     d1= mB[2][0]*mV[1] + mB[2][1]*mV[2] + mB[2][2]*mV[4] - mC[16];
2773     d2= mB[2][0]*mV[3] + mB[2][1]*mV[4] + mB[2][2]*mV[5] - mC[17];
2774
2775     //fC[15]+= mC[15]+ d0;
2776     //fC[16]+= mC[16]+ d1;
2777     //fC[17]+= mC[17]+ d2;
2778     fC[18]+= mC[18]+ d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
2779     fC[19]+= mC[19]+ d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
2780     fC[20]+= mC[20]+ d0*mB[2][0] + d1*mB[2][1] + d2*mB[2][2];
2781
2782     d0= mB[3][0]*mV[0] + mB[3][1]*mV[1] + mB[3][2]*mV[3] - mC[21];
2783     d1= mB[3][0]*mV[1] + mB[3][1]*mV[2] + mB[3][2]*mV[4] - mC[22];
2784     d2= mB[3][0]*mV[3] + mB[3][1]*mV[4] + mB[3][2]*mV[5] - mC[23];
2785
2786     //fC[21]+= mC[21] + d0;
2787     //fC[22]+= mC[22] + d1;
2788     //fC[23]+= mC[23] + d2;
2789     fC[24]+= mC[24] + d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
2790     fC[25]+= mC[25] + d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
2791     fC[26]+= mC[26] + d0*mB[2][0] + d1*mB[2][1] + d2*mB[2][2];
2792     fC[27]+= mC[27] + d0*mB[3][0] + d1*mB[3][1] + d2*mB[3][2];
2793   }
2794
2795 //  SetMassConstraint(0,0);
2796   SetNonlinearMassConstraint(0);
2797 }
2798
2799 Bool_t AliKFParticleBase::InvertSym3( const Double_t A[], Double_t Ai[] )
2800 {
2801   //* Invert symmetric matric stored in low-triagonal form 
2802
2803   bool ret = 0;
2804   double a0 = A[0], a1 = A[1], a2 = A[2], a3 = A[3];
2805
2806   Ai[0] = a2*A[5] - A[4]*A[4];
2807   Ai[1] = a3*A[4] - a1*A[5];
2808   Ai[3] = a1*A[4] - a2*a3;
2809   Double_t det = (a0*Ai[0] + a1*Ai[1] + a3*Ai[3]);
2810   if( TMath::Abs(det)>1.e-20 ) det = 1./det;    
2811   else{ 
2812     det = 0;
2813     ret = 1;
2814   }
2815   Ai[0] *= det;
2816   Ai[1] *= det;
2817   Ai[3] *= det;
2818   Ai[2] = ( a0*A[5] - a3*a3 )*det;
2819   Ai[4] = ( a1*a3 - a0*A[4] )*det;
2820   Ai[5] = ( a0*a2 - a1*a1 )*det;
2821   return ret;
2822 }
2823
2824 void AliKFParticleBase::MultQSQt( const Double_t Q[], const Double_t S[], Double_t SOut[] )
2825 {
2826   //* Matrix multiplication Q*S*Q^T, Q - square matrix, S - symmetric
2827
2828   const Int_t kN= 8;
2829   Double_t mA[kN*kN];
2830   
2831   for( Int_t i=0, ij=0; i<kN; i++ ){
2832     for( Int_t j=0; j<kN; j++, ++ij ){
2833       mA[ij] = 0 ;
2834       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];
2835     }
2836   }
2837     
2838   for( Int_t i=0; i<kN; i++ ){
2839     for( Int_t j=0; j<=i; j++ ){
2840       Int_t ij = ( j<=i ) ? i*(i+1)/2+j :j*(j+1)/2+i;
2841       SOut[ij] = 0 ;
2842       for( Int_t k=0; k<kN; k++ )  SOut[ij] += Q[ i*kN+k ] * mA[ k*kN+j ];
2843     }
2844   }
2845 }
2846
2847
2848 // 72-charachters line to define the printer border
2849 //3456789012345678901234567890123456789012345678901234567890123456789012
2850