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