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