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