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