]> git.uio.no Git - u/mrichter/AliRoot.git/blame - STEER/AliKFParticleBase.cxx
Initialisation of fractions corrected.
[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
5fc72f28 105Int_t AliKFParticleBase::GetMomentum( Double_t &p, Double_t &error ) const
f826d409 106{
107 //* Calculate particle momentum
108
109 Double_t x = fP[3];
110 Double_t y = fP[4];
111 Double_t z = fP[5];
55ac3e1e 112
f826d409 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);
55ac3e1e 118
5fc72f28 119 error = (x2*fC[9]+y2*fC[14]+z2*fC[20] + 2*(x*y*fC[13]+x*z*fC[18]+y*z*fC[19]) );
55ac3e1e 120 if( error>1.e-16 && p>1.e-4 ){
5fc72f28 121 error = TMath::Sqrt(error)/p;
f826d409 122 return 0;
123 }
55ac3e1e 124 error = 1.e8;
f826d409 125 return 1;
126}
127
5fc72f28 128Int_t AliKFParticleBase::GetPt( Double_t &pt, Double_t &error ) const
706952f5 129{
130 //* Calculate particle transverse momentum
131
132 Double_t px = fP[3];
133 Double_t py = fP[4];
134 Double_t px2 = px*px;
135 Double_t py2 = py*py;
136 Double_t pt2 = px2+py2;
5fc72f28 137 pt = TMath::Sqrt(pt2);
138 error = (px2*fC[9] + py2*fC[14] + 2*px*py*fC[13] );
139 if( error>0 && pt>1.e-4 ){
140 error = TMath::Sqrt(error)/pt;
706952f5 141 return 0;
142 }
5fc72f28 143 error = 1.e10;
706952f5 144 return 1;
145}
146
5fc72f28 147Int_t AliKFParticleBase::GetEta( Double_t &eta, Double_t &error ) const
706952f5 148{
149 //* Calculate particle pseudorapidity
150
151 Double_t px = fP[3];
152 Double_t py = fP[4];
153 Double_t pz = fP[5];
154 Double_t pt2 = px*px + py*py;
155 Double_t p2 = pt2 + pz*pz;
156 Double_t p = TMath::Sqrt(p2);
157 Double_t a = p + pz;
158 Double_t b = p - pz;
5fc72f28 159 eta = 1.e10;
706952f5 160 if( b > 1.e-8 ){
161 Double_t c = a/b;
a65041d0 162 if( c>1.e-8 ) eta = 0.5*TMath::Log(c);
706952f5 163 }
164 Double_t h3 = -px*pz;
a65041d0 165 Double_t h4 = -py*pz;
166 Double_t pt4 = pt2*pt2;
167 Double_t p2pt4 = p2*pt4;
168 error = (h3*h3*fC[9] + h4*h4*fC[14] + pt4*fC[20] + 2*( h3*(h4*fC[13] + fC[18]*pt2) + pt2*h4*fC[19] ) );
706952f5 169
a65041d0 170 if( error>0 && p2pt4>1.e-10 ){
171 error = TMath::Sqrt(error/p2pt4);
706952f5 172 return 0;
173 }
a65041d0 174
5fc72f28 175 error = 1.e10;
706952f5 176 return 1;
177}
178
5fc72f28 179Int_t AliKFParticleBase::GetPhi( Double_t &phi, Double_t &error ) const
706952f5 180{
181 //* Calculate particle polar angle
182
183 Double_t px = fP[3];
184 Double_t py = fP[4];
185 Double_t px2 = px*px;
186 Double_t py2 = py*py;
187 Double_t pt2 = px2 + py2;
5fc72f28 188 phi = TMath::ATan2(py,px);
189 error = (py2*fC[9] + px2*fC[14] - 2*px*py*fC[13] );
190 if( error>0 && pt2>1.e-4 ){
191 error = TMath::Sqrt(error)/pt2;
706952f5 192 return 0;
193 }
5fc72f28 194 error = 1.e10;
706952f5 195 return 1;
196}
197
5fc72f28 198Int_t AliKFParticleBase::GetR( Double_t &r, Double_t &error ) const
706952f5 199{
200 //* Calculate distance to the origin
201
202 Double_t x = fP[0];
203 Double_t y = fP[1];
204 Double_t x2 = x*x;
205 Double_t y2 = y*y;
5fc72f28 206 r = TMath::Sqrt(x2 + y2);
207 error = (x2*fC[0] + y2*fC[2] - 2*x*y*fC[1] );
208 if( error>0 && r>1.e-4 ){
209 error = TMath::Sqrt(error)/r;
706952f5 210 return 0;
211 }
5fc72f28 212 error = 1.e10;
706952f5 213 return 1;
214}
215
5fc72f28 216Int_t AliKFParticleBase::GetMass( Double_t &m, Double_t &error ) const
f826d409 217{
218 //* Calculate particle mass
219
220 // s = sigma^2 of m2/2
221
222 Double_t s = ( fP[3]*fP[3]*fC[9] + fP[4]*fP[4]*fC[14] + fP[5]*fP[5]*fC[20]
223 + fP[6]*fP[6]*fC[27]
224 +2*( + fP[3]*fP[4]*fC[13] + fP[5]*(fP[3]*fC[18] + fP[4]*fC[19])
225 - fP[6]*( fP[3]*fC[24] + fP[4]*fC[25] + fP[5]*fC[26] ) )
226 );
706952f5 227 Double_t m2 = TMath::Abs(fP[6]*fP[6] - fP[3]*fP[3] - fP[4]*fP[4] - fP[5]*fP[5]);
5fc72f28 228 m = TMath::Sqrt(m2);
229 if( m>1.e-10 ){
616ffc76 230 if( s>=0 ){
5fc72f28 231 error = TMath::Sqrt(s)/m;
616ffc76 232 return 0;
233 }
f826d409 234 }
5fc72f28 235 error = 1.e20;
f826d409 236 return 1;
237}
238
239
5fc72f28 240Int_t AliKFParticleBase::GetDecayLength( Double_t &l, Double_t &error ) const
f826d409 241{
242 //* Calculate particle decay length [cm]
243
244 Double_t x = fP[3];
245 Double_t y = fP[4];
246 Double_t z = fP[5];
247 Double_t t = fP[7];
248 Double_t x2 = x*x;
249 Double_t y2 = y*y;
250 Double_t z2 = z*z;
251 Double_t p2 = x2+y2+z2;
5fc72f28 252 l = t*TMath::Sqrt(p2);
f826d409 253 if( p2>1.e-4){
5fc72f28 254 error = p2*fC[35] + t*t/p2*(x2*fC[9]+y2*fC[14]+z2*fC[20]
f826d409 255 + 2*(x*y*fC[13]+x*z*fC[18]+y*z*fC[19]) )
256 + 2*t*(x*fC[31]+y*fC[32]+z*fC[33]);
5fc72f28 257 error = TMath::Sqrt(TMath::Abs(error));
f826d409 258 return 0;
259 }
5fc72f28 260 error = 1.e20;
f826d409 261 return 1;
262}
263
446ce366 264Int_t AliKFParticleBase::GetDecayLengthXY( Double_t &l, Double_t &error ) const
265{
266 //* Calculate particle decay length in XY projection [cm]
267
268 Double_t x = fP[3];
269 Double_t y = fP[4];
270 Double_t t = fP[7];
271 Double_t x2 = x*x;
272 Double_t y2 = y*y;
273 Double_t pt2 = x2+y2;
274 l = t*TMath::Sqrt(pt2);
275 if( pt2>1.e-4){
276 error = pt2*fC[35] + t*t/pt2*(x2*fC[9]+y2*fC[14] + 2*x*y*fC[13] )
277 + 2*t*(x*fC[31]+y*fC[32]);
278 error = TMath::Sqrt(TMath::Abs(error));
279 return 0;
280 }
281 error = 1.e20;
282 return 1;
283}
284
285
5fc72f28 286Int_t AliKFParticleBase::GetLifeTime( Double_t &tauC, Double_t &error ) const
f826d409 287{
288 //* Calculate particle decay time [s]
289
290 Double_t m, dm;
291 GetMass( m, dm );
292 Double_t cTM = (-fP[3]*fC[31] - fP[4]*fC[32] - fP[5]*fC[33] + fP[6]*fC[34]);
5fc72f28 293 tauC = fP[7]*m;
294 error = m*m*fC[35] + 2*fP[7]*cTM + fP[7]*fP[7]*dm*dm;
295 if( error > 0 ){
296 error = TMath::Sqrt( error );
f826d409 297 return 0;
298 }
5fc72f28 299 error = 1.e20;
f826d409 300 return 1;
301}
302
303
304void AliKFParticleBase::operator +=( const AliKFParticleBase &Daughter )
305{
306 //* Add daughter via operator+=
307
308 AddDaughter( Daughter );
309}
310
706952f5 311Double_t AliKFParticleBase::GetSCorrection( const Double_t Part[], const Double_t XYZ[] )
312{
313 //* Get big enough correction for S error to let the particle Part be fitted to XYZ point
314
315 Double_t d[3] = { XYZ[0]-Part[0], XYZ[1]-Part[1], XYZ[2]-Part[2] };
316 Double_t p2 = Part[3]*Part[3]+Part[4]*Part[4]+Part[5]*Part[5];
a65041d0 317 Double_t sigmaS = (p2>1.e-4) ? ( 10.1+3.*TMath::Sqrt( d[0]*d[0]+d[1]*d[1]+d[2]*d[2]) )/TMath::Sqrt(p2) : 1.;
706952f5 318 return sigmaS;
319}
616ffc76 320
321void AliKFParticleBase::GetMeasurement( const Double_t XYZ[], Double_t m[], Double_t V[] ) const
322{
e7b09c95 323 //* Get additional covariances V used during measurement
324
616ffc76 325 Double_t b[3];
326 GetFieldValue( XYZ, b );
327 const Double_t kCLight = 0.000299792458;
328 b[0]*=kCLight; b[1]*=kCLight; b[2]*=kCLight;
329
330 Transport( GetDStoPoint(XYZ), m, V );
331
706952f5 332 Double_t sigmaS = GetSCorrection( m, XYZ );
616ffc76 333
334 Double_t h[6];
335
336 h[0] = m[3]*sigmaS;
337 h[1] = m[4]*sigmaS;
338 h[2] = m[5]*sigmaS;
339 h[3] = ( h[1]*b[2]-h[2]*b[1] )*GetQ();
340 h[4] = ( h[2]*b[0]-h[0]*b[2] )*GetQ();
341 h[5] = ( h[0]*b[1]-h[1]*b[0] )*GetQ();
342
616ffc76 343 V[ 0]+= h[0]*h[0];
344 V[ 1]+= h[1]*h[0];
345 V[ 2]+= h[1]*h[1];
346 V[ 3]+= h[2]*h[0];
347 V[ 4]+= h[2]*h[1];
348 V[ 5]+= h[2]*h[2];
349
350 V[ 6]+= h[3]*h[0];
351 V[ 7]+= h[3]*h[1];
352 V[ 8]+= h[3]*h[2];
353 V[ 9]+= h[3]*h[3];
354
355 V[10]+= h[4]*h[0];
356 V[11]+= h[4]*h[1];
357 V[12]+= h[4]*h[2];
358 V[13]+= h[4]*h[3];
359 V[14]+= h[4]*h[4];
360
361 V[15]+= h[5]*h[0];
362 V[16]+= h[5]*h[1];
363 V[17]+= h[5]*h[2];
364 V[18]+= h[5]*h[3];
365 V[19]+= h[5]*h[4];
366 V[20]+= h[5]*h[5];
367}
368
616ffc76 369
f826d409 370void AliKFParticleBase::AddDaughter( const AliKFParticleBase &Daughter )
371{
372 //* Add daughter
373
616ffc76 374 if( fNDF<-1 ){ // first daughter -> just copy
375 fNDF = -1;
376 fQ = Daughter.GetQ();
377 for( Int_t i=0; i<7; i++) fP[i] = Daughter.fP[i];
378 for( Int_t i=0; i<28; i++) fC[i] = Daughter.fC[i];
379 fSFromDecay = 0;
380 return;
381 }
382
f826d409 383 TransportToDecayVertex();
384
385 Double_t b[3];
386 Int_t maxIter = 1;
387
388 if( !fIsLinearized ){
389 if( fNDF==-1 ){
390 Double_t ds, ds1;
616ffc76 391 GetDStoParticle(Daughter, ds, ds1);
f826d409 392 TransportToDS( ds );
616ffc76 393 Double_t m[8];
394 Double_t mCd[36];
395 Daughter.Transport( ds1, m, mCd );
396 fVtxGuess[0] = .5*( fP[0] + m[0] );
397 fVtxGuess[1] = .5*( fP[1] + m[1] );
398 fVtxGuess[2] = .5*( fP[2] + m[2] );
399 } else {
400 fVtxGuess[0] = fP[0];
401 fVtxGuess[1] = fP[1];
402 fVtxGuess[2] = fP[2];
f826d409 403 }
f826d409 404 maxIter = 3;
405 }
406
407 for( Int_t iter=0; iter<maxIter; iter++ ){
408
409 {
410 GetFieldValue( fVtxGuess, b );
411 const Double_t kCLight = 0.000299792458;
412 b[0]*=kCLight; b[1]*=kCLight; b[2]*=kCLight;
413 }
f826d409 414
616ffc76 415 Double_t *ffP = fP, *ffC = fC, tmpP[8], tmpC[36];
de0d0ceb 416 if( fNDF==-1 ){
616ffc76 417 GetMeasurement( fVtxGuess, tmpP, tmpC );
418 ffP = tmpP;
419 ffC = tmpC;
f826d409 420 }
421
e7b09c95 422 Double_t m[8], mV[36];
de0d0ceb 423
e7b09c95 424 if( Daughter.fC[35]>0 ){
425 Daughter.GetMeasurement( fVtxGuess, m, mV );
426 } else {
427 for( Int_t i=0; i<8; i++ ) m[i] = Daughter.fP[i];
428 for( Int_t i=0; i<36; i++ ) mV[i] = Daughter.fC[i];
429 }
de0d0ceb 430
f826d409 431 //*
432
433 Double_t mS[6];
434 {
616ffc76 435 Double_t mSi[6] = { ffC[0]+mV[0],
436 ffC[1]+mV[1], ffC[2]+mV[2],
437 ffC[3]+mV[3], ffC[4]+mV[4], ffC[5]+mV[5] };
de0d0ceb 438
f826d409 439 mS[0] = mSi[2]*mSi[5] - mSi[4]*mSi[4];
440 mS[1] = mSi[3]*mSi[4] - mSi[1]*mSi[5];
441 mS[2] = mSi[0]*mSi[5] - mSi[3]*mSi[3];
442 mS[3] = mSi[1]*mSi[4] - mSi[2]*mSi[3];
443 mS[4] = mSi[1]*mSi[3] - mSi[0]*mSi[4];
444 mS[5] = mSi[0]*mSi[2] - mSi[1]*mSi[1];
445
706952f5 446 Double_t s = ( mSi[0]*mS[0] + mSi[1]*mS[1] + mSi[3]*mS[3] );
447
f826d409 448 s = ( s > 1.E-20 ) ?1./s :0;
449 mS[0]*=s;
450 mS[1]*=s;
451 mS[2]*=s;
452 mS[3]*=s;
453 mS[4]*=s;
454 mS[5]*=s;
455 }
456
457 //* Residual (measured - estimated)
458
616ffc76 459 Double_t zeta[3] = { m[0]-ffP[0], m[1]-ffP[1], m[2]-ffP[2] };
de0d0ceb 460
f826d409 461
462 //* CHt = CH' - D'
463
464 Double_t mCHt0[7], mCHt1[7], mCHt2[7];
465
616ffc76 466 mCHt0[0]=ffC[ 0] ; mCHt1[0]=ffC[ 1] ; mCHt2[0]=ffC[ 3] ;
467 mCHt0[1]=ffC[ 1] ; mCHt1[1]=ffC[ 2] ; mCHt2[1]=ffC[ 4] ;
468 mCHt0[2]=ffC[ 3] ; mCHt1[2]=ffC[ 4] ; mCHt2[2]=ffC[ 5] ;
469 mCHt0[3]=ffC[ 6]-mV[ 6]; mCHt1[3]=ffC[ 7]-mV[ 7]; mCHt2[3]=ffC[ 8]-mV[ 8];
470 mCHt0[4]=ffC[10]-mV[10]; mCHt1[4]=ffC[11]-mV[11]; mCHt2[4]=ffC[12]-mV[12];
471 mCHt0[5]=ffC[15]-mV[15]; mCHt1[5]=ffC[16]-mV[16]; mCHt2[5]=ffC[17]-mV[17];
472 mCHt0[6]=ffC[21]-mV[21]; mCHt1[6]=ffC[22]-mV[22]; mCHt2[6]=ffC[23]-mV[23];
f826d409 473
474 //* Kalman gain K = mCH'*S
475
476 Double_t k0[7], k1[7], k2[7];
477
478 for(Int_t i=0;i<7;++i){
479 k0[i] = mCHt0[i]*mS[0] + mCHt1[i]*mS[1] + mCHt2[i]*mS[3];
480 k1[i] = mCHt0[i]*mS[1] + mCHt1[i]*mS[2] + mCHt2[i]*mS[4];
481 k2[i] = mCHt0[i]*mS[3] + mCHt1[i]*mS[4] + mCHt2[i]*mS[5];
482 }
483
4bbc290d 484 //* New estimation of the vertex position
f826d409 485
486 if( iter<maxIter-1 ){
487 for(Int_t i=0; i<3; ++i)
616ffc76 488 fVtxGuess[i]= ffP[i] + k0[i]*zeta[0]+k1[i]*zeta[1]+k2[i]*zeta[2];
f826d409 489 continue;
490 }
491
492 // last itearation -> update the particle
493
494 //* Add the daughter momentum to the particle momentum
495
616ffc76 496 ffP[ 3] += m[ 3];
497 ffP[ 4] += m[ 4];
498 ffP[ 5] += m[ 5];
499 ffP[ 6] += m[ 6];
f826d409 500
616ffc76 501 ffC[ 9] += mV[ 9];
502 ffC[13] += mV[13];
503 ffC[14] += mV[14];
504 ffC[18] += mV[18];
505 ffC[19] += mV[19];
506 ffC[20] += mV[20];
507 ffC[24] += mV[24];
508 ffC[25] += mV[25];
509 ffC[26] += mV[26];
510 ffC[27] += mV[27];
f826d409 511
de0d0ceb 512
513 //* New estimation of the vertex position r += K*zeta
f826d409 514
515 for(Int_t i=0;i<7;++i)
616ffc76 516 fP[i] = ffP[i] + k0[i]*zeta[0] + k1[i]*zeta[1] + k2[i]*zeta[2];
f826d409 517
518 //* New covariance matrix C -= K*(mCH')'
de0d0ceb 519
f826d409 520 for(Int_t i=0, k=0;i<7;++i){
de0d0ceb 521 for(Int_t j=0;j<=i;++j,++k){
616ffc76 522 fC[k] = ffC[k] - (k0[i]*mCHt0[j] + k1[i]*mCHt1[j] + k2[i]*mCHt2[j] );
de0d0ceb 523 }
f826d409 524 }
de0d0ceb 525
f826d409 526 //* Calculate Chi^2
527
528 fNDF += 2;
529 fQ += Daughter.GetQ();
530 fSFromDecay = 0;
531 fChi2 += (mS[0]*zeta[0] + mS[1]*zeta[1] + mS[3]*zeta[2])*zeta[0]
532 + (mS[1]*zeta[0] + mS[2]*zeta[1] + mS[4]*zeta[2])*zeta[1]
706952f5 533 + (mS[3]*zeta[0] + mS[4]*zeta[1] + mS[5]*zeta[2])*zeta[2];
de0d0ceb 534
f826d409 535 }
536}
537
f826d409 538void AliKFParticleBase::SetProductionVertex( const AliKFParticleBase &Vtx )
539{
706952f5 540 //* Set production vertex for the particle, when the particle was not used in the vertex fit
f826d409 541
542 const Double_t *m = Vtx.fP, *mV = Vtx.fC;
543
e7b09c95 544 Bool_t noS = ( fC[35]<=0 ); // no decay length allowed
545
546 if( noS ){
547 TransportToDecayVertex();
548 fP[7] = 0;
55ac3e1e 549 fC[28] = fC[29] = fC[30] = fC[31] = fC[32] = fC[33] = fC[34] = fC[35] = 0;
e7b09c95 550 } else {
55ac3e1e 551 TransportToDS( GetDStoPoint( m ) );
e7b09c95 552 fP[7] = -fSFromDecay;
55ac3e1e 553 fC[28] = fC[29] = fC[30] = fC[31] = fC[32] = fC[33] = fC[34] = 0;
554 fC[35] = 1000.;
555
e7b09c95 556 Convert(1);
557 }
558
f826d409 559 Double_t mAi[6];
706952f5 560
446ce366 561 InvertSym3( fC, mAi );
f826d409 562
563 Double_t mB[5][3];
564
565 mB[0][0] = fC[ 6]*mAi[0] + fC[ 7]*mAi[1] + fC[ 8]*mAi[3];
566 mB[0][1] = fC[ 6]*mAi[1] + fC[ 7]*mAi[2] + fC[ 8]*mAi[4];
567 mB[0][2] = fC[ 6]*mAi[3] + fC[ 7]*mAi[4] + fC[ 8]*mAi[5];
568
569 mB[1][0] = fC[10]*mAi[0] + fC[11]*mAi[1] + fC[12]*mAi[3];
570 mB[1][1] = fC[10]*mAi[1] + fC[11]*mAi[2] + fC[12]*mAi[4];
571 mB[1][2] = fC[10]*mAi[3] + fC[11]*mAi[4] + fC[12]*mAi[5];
572
573 mB[2][0] = fC[15]*mAi[0] + fC[16]*mAi[1] + fC[17]*mAi[3];
574 mB[2][1] = fC[15]*mAi[1] + fC[16]*mAi[2] + fC[17]*mAi[4];
575 mB[2][2] = fC[15]*mAi[3] + fC[16]*mAi[4] + fC[17]*mAi[5];
576
577 mB[3][0] = fC[21]*mAi[0] + fC[22]*mAi[1] + fC[23]*mAi[3];
578 mB[3][1] = fC[21]*mAi[1] + fC[22]*mAi[2] + fC[23]*mAi[4];
579 mB[3][2] = fC[21]*mAi[3] + fC[22]*mAi[4] + fC[23]*mAi[5];
580
581 mB[4][0] = fC[28]*mAi[0] + fC[29]*mAi[1] + fC[30]*mAi[3];
582 mB[4][1] = fC[28]*mAi[1] + fC[29]*mAi[2] + fC[30]*mAi[4];
583 mB[4][2] = fC[28]*mAi[3] + fC[29]*mAi[4] + fC[30]*mAi[5];
584
585 Double_t z[3] = { m[0]-fP[0], m[1]-fP[1], m[2]-fP[2] };
586
587 {
a65041d0 588 Double_t mAVi[6] = { fC[0]-mV[0], fC[1]-mV[1], fC[2]-mV[2],
f826d409 589 fC[3]-mV[3], fC[4]-mV[4], fC[5]-mV[5] };
706952f5 590
a65041d0 591 if( !InvertSym3( mAVi, mAVi ) ){
706952f5 592
593 Double_t dChi2 = ( +(mAVi[0]*z[0] + mAVi[1]*z[1] + mAVi[3]*z[2])*z[0]
594 +(mAVi[1]*z[0] + mAVi[2]*z[1] + mAVi[4]*z[2])*z[1]
a65041d0 595 +(mAVi[3]*z[0] + mAVi[4]*z[1] + mAVi[5]*z[2])*z[2] );
706952f5 596
597 // Take Abs(dChi2) here. Negative value of 'det' or 'dChi2' shows that the particle
598 // was not used in the production vertex fit
599
600 fChi2+= TMath::Abs( dChi2 );
601 }
f826d409 602 fNDF += 2;
f826d409 603 }
604
605 fP[0] = m[0];
606 fP[1] = m[1];
607 fP[2] = m[2];
608 fP[3]+= mB[0][0]*z[0] + mB[0][1]*z[1] + mB[0][2]*z[2];
609 fP[4]+= mB[1][0]*z[0] + mB[1][1]*z[1] + mB[1][2]*z[2];
610 fP[5]+= mB[2][0]*z[0] + mB[2][1]*z[1] + mB[2][2]*z[2];
611 fP[6]+= mB[3][0]*z[0] + mB[3][1]*z[1] + mB[3][2]*z[2];
612 fP[7]+= mB[4][0]*z[0] + mB[4][1]*z[1] + mB[4][2]*z[2];
613
614 Double_t d0, d1, d2;
615
616 fC[0] = mV[0];
617 fC[1] = mV[1];
618 fC[2] = mV[2];
619 fC[3] = mV[3];
620 fC[4] = mV[4];
621 fC[5] = mV[5];
622
623 d0= mB[0][0]*mV[0] + mB[0][1]*mV[1] + mB[0][2]*mV[3] - fC[ 6];
624 d1= mB[0][0]*mV[1] + mB[0][1]*mV[2] + mB[0][2]*mV[4] - fC[ 7];
625 d2= mB[0][0]*mV[3] + mB[0][1]*mV[4] + mB[0][2]*mV[5] - fC[ 8];
626
627 fC[ 6]+= d0;
628 fC[ 7]+= d1;
629 fC[ 8]+= d2;
630 fC[ 9]+= d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
631
632 d0= mB[1][0]*mV[0] + mB[1][1]*mV[1] + mB[1][2]*mV[3] - fC[10];
633 d1= mB[1][0]*mV[1] + mB[1][1]*mV[2] + mB[1][2]*mV[4] - fC[11];
634 d2= mB[1][0]*mV[3] + mB[1][1]*mV[4] + mB[1][2]*mV[5] - fC[12];
635
636 fC[10]+= d0;
637 fC[11]+= d1;
638 fC[12]+= d2;
639 fC[13]+= d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
640 fC[14]+= d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
641
642 d0= mB[2][0]*mV[0] + mB[2][1]*mV[1] + mB[2][2]*mV[3] - fC[15];
643 d1= mB[2][0]*mV[1] + mB[2][1]*mV[2] + mB[2][2]*mV[4] - fC[16];
644 d2= mB[2][0]*mV[3] + mB[2][1]*mV[4] + mB[2][2]*mV[5] - fC[17];
645
646 fC[15]+= d0;
647 fC[16]+= d1;
648 fC[17]+= d2;
649 fC[18]+= d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
650 fC[19]+= d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
651 fC[20]+= d0*mB[2][0] + d1*mB[2][1] + d2*mB[2][2];
652
653 d0= mB[3][0]*mV[0] + mB[3][1]*mV[1] + mB[3][2]*mV[3] - fC[21];
654 d1= mB[3][0]*mV[1] + mB[3][1]*mV[2] + mB[3][2]*mV[4] - fC[22];
655 d2= mB[3][0]*mV[3] + mB[3][1]*mV[4] + mB[3][2]*mV[5] - fC[23];
656
657 fC[21]+= d0;
658 fC[22]+= d1;
659 fC[23]+= d2;
660 fC[24]+= d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
661 fC[25]+= d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
662 fC[26]+= d0*mB[2][0] + d1*mB[2][1] + d2*mB[2][2];
663 fC[27]+= d0*mB[3][0] + d1*mB[3][1] + d2*mB[3][2];
664
665 d0= mB[4][0]*mV[0] + mB[4][1]*mV[1] + mB[4][2]*mV[3] - fC[28];
666 d1= mB[4][0]*mV[1] + mB[4][1]*mV[2] + mB[4][2]*mV[4] - fC[29];
667 d2= mB[4][0]*mV[3] + mB[4][1]*mV[4] + mB[4][2]*mV[5] - fC[30];
668
669 fC[28]+= d0;
670 fC[29]+= d1;
671 fC[30]+= d2;
672 fC[31]+= d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
673 fC[32]+= d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
674 fC[33]+= d0*mB[2][0] + d1*mB[2][1] + d2*mB[2][2];
675 fC[34]+= d0*mB[3][0] + d1*mB[3][1] + d2*mB[3][2];
676 fC[35]+= d0*mB[4][0] + d1*mB[4][1] + d2*mB[4][2];
677
e7b09c95 678 if( noS ){
679 fP[7] = 0;
55ac3e1e 680 fC[28] = fC[29] = fC[30] = fC[31] = fC[32] = fC[33] = fC[34] = fC[35] = 0;
e7b09c95 681 } else {
682 TransportToDS( fP[7] );
683 Convert(0);
684 }
685
f826d409 686 fSFromDecay = 0;
f826d409 687}
688
689
690
e7b09c95 691void AliKFParticleBase::SetMassConstraint( Double_t Mass, Double_t SigmaMass )
f826d409 692{
706952f5 693 //* Set hard( SigmaMass=0 ) or soft (SigmaMass>0) mass constraint
694
695 Double_t m2 = Mass*Mass; // measurement, weighted by Mass
696 Double_t s2 = m2*SigmaMass*SigmaMass; // sigma^2
697
698 Double_t p2 = fP[3]*fP[3] + fP[4]*fP[4] + fP[5]*fP[5];
699 Double_t e0 = TMath::Sqrt(m2+p2);
f826d409 700
701 Double_t mH[8];
702 mH[0] = mH[1] = mH[2] = 0.;
703 mH[3] = -2*fP[3];
704 mH[4] = -2*fP[4];
705 mH[5] = -2*fP[5];
706952f5 706 mH[6] = 2*fP[6];//e0;
f826d409 707 mH[7] = 0;
f826d409 708
706952f5 709 Double_t zeta = e0*e0 - e0*fP[6];
710 zeta = m2 - (fP[6]*fP[6]-p2);
f826d409 711
706952f5 712 Double_t mCHt[8], s2_est=0;
713 for( Int_t i=0; i<8; ++i ){
f826d409 714 mCHt[i] = 0.0;
715 for (Int_t j=0;j<8;++j) mCHt[i] += Cij(i,j)*mH[j];
706952f5 716 s2_est += mH[i]*mCHt[i];
f826d409 717 }
718
706952f5 719 if( s2_est<1.e-20 ) return; // calculated mass error is already 0,
720 // the particle can not be constrained on mass
721
722 Double_t w2 = 1./( s2 + s2_est );
723 fChi2 += zeta*zeta*w2;
f826d409 724 fNDF += 1;
725 for( Int_t i=0, ii=0; i<8; ++i ){
706952f5 726 Double_t ki = mCHt[i]*w2;
f826d409 727 fP[i]+= ki*zeta;
728 for(Int_t j=0;j<=i;++j) fC[ii++] -= ki*mCHt[j];
729 }
730}
731
732
e7b09c95 733void AliKFParticleBase::SetNoDecayLength()
734{
735 //* Set no decay length for resonances
736
737 TransportToDecayVertex();
738
739 Double_t h[8];
740 h[0] = h[1] = h[2] = h[3] = h[4] = h[5] = h[6] = 0;
741 h[7] = 1;
742
743 Double_t zeta = 0 - fP[7];
744 for(Int_t i=0;i<8;++i) zeta -= h[i]*(fP[i]-fP[i]);
745
746 Double_t s = fC[35];
747 if( s>1.e-20 ){
748 s = 1./s;
749 fChi2 += zeta*zeta*s;
750 fNDF += 1;
751 for( Int_t i=0, ii=0; i<7; ++i ){
752 Double_t ki = fC[28+i]*s;
753 fP[i]+= ki*zeta;
754 for(Int_t j=0;j<=i;++j) fC[ii++] -= ki*fC[28+j];
755 }
756 }
757 fP[7] = 0;
55ac3e1e 758 fC[28] = fC[29] = fC[30] = fC[31] = fC[32] = fC[33] = fC[34] = fC[35] = 0;
e7b09c95 759}
760
761
f826d409 762void AliKFParticleBase::Construct( const AliKFParticleBase* vDaughters[], Int_t NDaughters,
706952f5 763 const AliKFParticleBase *Parent, Double_t Mass, Bool_t IsConstrained )
f826d409 764{
765 //* Full reconstruction in one go
766
767 Int_t maxIter = 1;
768 bool wasLinearized = fIsLinearized;
706952f5 769 if( !fIsLinearized || IsConstrained ){
616ffc76 770 //fVtxGuess[0] = fVtxGuess[1] = fVtxGuess[2] = 0; //!!!!
771 fVtxGuess[0] = GetX();
772 fVtxGuess[1] = GetY();
773 fVtxGuess[2] = GetZ();
f826d409 774 fIsLinearized = 1;
775 maxIter = 3;
776 }
777
706952f5 778 Double_t constraintC[6];
779
780 if( IsConstrained ){
781 for(Int_t i=0;i<6;++i) constraintC[i]=fC[i];
782 } else {
783 for(Int_t i=0;i<6;++i) constraintC[i]=0.;
784 constraintC[0] = constraintC[2] = constraintC[5] = 100.;
785 }
786
787
f826d409 788 for( Int_t iter=0; iter<maxIter; iter++ ){
f826d409 789 fAtProductionVertex = 0;
790 fSFromDecay = 0;
791 fP[0] = fVtxGuess[0];
792 fP[1] = fVtxGuess[1];
793 fP[2] = fVtxGuess[2];
794 fP[3] = 0;
795 fP[4] = 0;
796 fP[5] = 0;
797 fP[6] = 0;
798 fP[7] = 0;
706952f5 799
800 for(Int_t i=0;i<6; ++i) fC[i]=constraintC[i];
801 for(Int_t i=6;i<36;++i) fC[i]=0.;
f826d409 802 fC[35] = 1.;
803
706952f5 804 fNDF = IsConstrained ?0 :-3;
f826d409 805 fChi2 = 0.;
806 fQ = 0;
807
4bbc290d 808 for( Int_t itr =0; itr<NDaughters; itr++ ){
809 AddDaughter( *vDaughters[itr] );
810 }
f826d409 811 if( iter<maxIter-1){
812 for( Int_t i=0; i<3; i++ ) fVtxGuess[i] = fP[i];
813 }
814 }
815 fIsLinearized = wasLinearized;
816
817 if( Mass>=0 ) SetMassConstraint( Mass );
818 if( Parent ) SetProductionVertex( *Parent );
819}
820
821
822void AliKFParticleBase::Convert( bool ToProduction )
823{
824 //* Tricky function - convert the particle error along its trajectory to
825 //* the value which corresponds to its production/decay vertex
826 //* It is done by combination of the error of decay length with the position errors
827
828 Double_t fld[3];
829 {
830 GetFieldValue( fP, fld );
831 const Double_t kCLight = fQ*0.000299792458;
832 fld[0]*=kCLight; fld[1]*=kCLight; fld[2]*=kCLight;
833 }
834
835 Double_t h[6];
836
837 h[0] = fP[3];
838 h[1] = fP[4];
839 h[2] = fP[5];
840 if( ToProduction ){ h[0]=-h[0]; h[1]=-h[1]; h[2]=-h[2]; }
841 h[3] = h[1]*fld[2]-h[2]*fld[1];
842 h[4] = h[2]*fld[0]-h[0]*fld[2];
843 h[5] = h[0]*fld[1]-h[1]*fld[0];
844
845 Double_t c;
846
847 c = fC[28]+h[0]*fC[35];
848 fC[ 0]+= h[0]*(c+fC[28]);
849 fC[28] = c;
850
851 fC[ 1]+= h[1]*fC[28] + h[0]*fC[29];
852 c = fC[29]+h[1]*fC[35];
853 fC[ 2]+= h[1]*(c+fC[29]);
854 fC[29] = c;
855
856 fC[ 3]+= h[2]*fC[28] + h[0]*fC[30];
857 fC[ 4]+= h[2]*fC[29] + h[1]*fC[30];
858 c = fC[30]+h[2]*fC[35];
859 fC[ 5]+= h[2]*(c+fC[30]);
860 fC[30] = c;
861
862 fC[ 6]+= h[3]*fC[28] + h[0]*fC[31];
863 fC[ 7]+= h[3]*fC[29] + h[1]*fC[31];
864 fC[ 8]+= h[3]*fC[30] + h[2]*fC[31];
865 c = fC[31]+h[3]*fC[35];
866 fC[ 9]+= h[3]*(c+fC[31]);
867 fC[31] = c;
868
869 fC[10]+= h[4]*fC[28] + h[0]*fC[32];
870 fC[11]+= h[4]*fC[29] + h[1]*fC[32];
871 fC[12]+= h[4]*fC[30] + h[2]*fC[32];
872 fC[13]+= h[4]*fC[31] + h[3]*fC[32];
873 c = fC[32]+h[4]*fC[35];
874 fC[14]+= h[4]*(c+fC[32]);
875 fC[32] = c;
876
877 fC[15]+= h[5]*fC[28] + h[0]*fC[33];
878 fC[16]+= h[5]*fC[29] + h[1]*fC[33];
879 fC[17]+= h[5]*fC[30] + h[2]*fC[33];
880 fC[18]+= h[5]*fC[31] + h[3]*fC[33];
881 fC[19]+= h[5]*fC[32] + h[4]*fC[33];
882 c = fC[33]+h[5]*fC[35];
883 fC[20]+= h[5]*(c+fC[33]);
884 fC[33] = c;
885
886 fC[21]+= h[0]*fC[34];
887 fC[22]+= h[1]*fC[34];
888 fC[23]+= h[2]*fC[34];
889 fC[24]+= h[3]*fC[34];
890 fC[25]+= h[4]*fC[34];
891 fC[26]+= h[5]*fC[34];
892}
893
894
895void AliKFParticleBase::TransportToDecayVertex()
896{
897 //* Transport the particle to its decay vertex
898
899 if( fSFromDecay != 0 ) TransportToDS( -fSFromDecay );
900 if( fAtProductionVertex ) Convert(0);
901 fAtProductionVertex = 0;
902}
903
904void AliKFParticleBase::TransportToProductionVertex()
905{
906 //* Transport the particle to its production vertex
616ffc76 907
f826d409 908 if( fSFromDecay != -fP[7] ) TransportToDS( -fSFromDecay-fP[7] );
909 if( !fAtProductionVertex ) Convert( 1 );
910 fAtProductionVertex = 1;
911}
912
913
914void AliKFParticleBase::TransportToDS( Double_t dS )
915{
916 //* Transport the particle on dS parameter (SignedPath/Momentum)
616ffc76 917
f826d409 918 Transport( dS, fP, fC );
919 fSFromDecay+= dS;
920}
921
922
923Double_t AliKFParticleBase::GetDStoPointLine( const Double_t xyz[] ) const
924{
925 //* Get dS to a certain space point without field
926
927 Double_t p2 = fP[3]*fP[3] + fP[4]*fP[4] + fP[5]*fP[5];
928 if( p2<1.e-4 ) p2 = 1;
929 return ( fP[3]*(xyz[0]-fP[0]) + fP[4]*(xyz[1]-fP[1]) + fP[5]*(xyz[2]-fP[2]) )/p2;
930}
931
932
933Double_t AliKFParticleBase::GetDStoPointBz( Double_t B, const Double_t xyz[] )
934 const
935{
616ffc76 936
f826d409 937 //* Get dS to a certain space point for Bz field
f826d409 938 const Double_t kCLight = 0.000299792458;
939 Double_t bq = B*fQ*kCLight;
940 Double_t pt2 = fP[3]*fP[3] + fP[4]*fP[4];
941 if( pt2<1.e-4 ) return 0;
942 Double_t dx = xyz[0] - fP[0];
943 Double_t dy = xyz[1] - fP[1];
944 Double_t a = dx*fP[3]+dy*fP[4];
706952f5 945 Double_t dS;
946
947 if( TMath::Abs(bq)<1.e-8 ) dS = a/pt2;
948 else dS = TMath::ATan2( bq*a, pt2 + bq*(dy*fP[3] -dx*fP[4]) )/bq;
949
950 if(0){
951
952 Double_t px = fP[3];
953 Double_t py = fP[4];
954 Double_t pz = fP[5];
955 Double_t ss[2], g[2][5];
956
957 ss[0] = dS;
958 ss[1] = -dS;
959 for( Int_t i=0; i<2; i++){
960 Double_t bs = bq*ss[i];
961 Double_t c = TMath::Cos(bs), s = TMath::Sin(bs);
962 Double_t cB,sB;
963 if( TMath::Abs(bq)>1.e-8){
964 cB= (1-c)/bq;
965 sB= s/bq;
966 }else{
bfd20868 967 const Double_t kOvSqr6 = 1./TMath::Sqrt(6.);
968 sB = (1.-bs*kOvSqr6)*(1.+bs*kOvSqr6)*ss[i];
706952f5 969 cB = .5*sB*bs;
970 }
971 g[i][0] = fP[0] + sB*px + cB*py;
972 g[i][1] = fP[1] - cB*px + sB*py;
973 g[i][2] = fP[2] + ss[i]*pz;
974 g[i][3] = + c*px + s*py;
975 g[i][4] = - s*px + c*py;
976 }
977
978 Int_t i=0;
979
980 Double_t dMin = 1.e10;
981 for( Int_t j=0; j<2; j++){
982 Double_t xx = g[j][0]-xyz[0];
983 Double_t yy = g[j][1]-xyz[1];
984 Double_t zz = g[j][2]-xyz[2];
985 Double_t d = xx*xx + yy*yy + zz*zz;
986 if( d<dMin ){
987 dMin = d;
988 i = j;
989 }
990 }
991
992 dS = ss[i];
993
994 Double_t x= g[i][0], y= g[i][1], z= g[i][2], ppx= g[i][3], ppy= g[i][4];
5fc72f28 995 Double_t ddx = x-xyz[0];
996 Double_t ddy = y-xyz[1];
997 Double_t ddz = z-xyz[2];
998 Double_t c = ddx*ppx + ddy*ppy + ddz*pz ;
706952f5 999 Double_t pp2 = ppx*ppx + ppy*ppy + pz*pz;
1000 if( TMath::Abs(pp2)>1.e-8 ){
1001 dS+=c/pp2;
1002 }
1003 }
1004 return dS;
f826d409 1005}
1006
1007
1008void AliKFParticleBase::GetDStoParticleBz( Double_t B, const AliKFParticleBase &p,
616ffc76 1009 Double_t &DS, Double_t &DS1 )
f826d409 1010 const
1011{
1012 //* Get dS to another particle for Bz field
f826d409 1013 Double_t px = fP[3];
1014 Double_t py = fP[4];
1015 Double_t pz = fP[5];
1016
1017 Double_t px1 = p.fP[3];
1018 Double_t py1 = p.fP[4];
1019 Double_t pz1 = p.fP[5];
1020
1021 const Double_t kCLight = 0.000299792458;
1022
1023 Double_t bq = B*fQ*kCLight;
1024 Double_t bq1 = B*p.fQ*kCLight;
f826d409 1025 Double_t s=0, ds=0, s1=0, ds1=0;
1026
1027 if( TMath::Abs(bq)>1.e-8 || TMath::Abs(bq1)>1.e-8 ){
4bbc290d 1028
f826d409 1029 Double_t dx = (p.fP[0] - fP[0]);
1030 Double_t dy = (p.fP[1] - fP[1]);
1031 Double_t d2 = (dx*dx+dy*dy);
1032
1033 Double_t p2 = (px *px + py *py);
1034 Double_t p21 = (px1*px1 + py1*py1);
1035
1036 Double_t a = (px*py1 - py*px1);
1037 Double_t b = (px*px1 + py*py1);
1038
1039 Double_t ldx = bq*bq1*dx - bq1*py + bq*py1 ;
1040 Double_t ldy = bq*bq1*dy + bq1*px - bq*px1 ;
1041 Double_t l2 = ldx*ldx + ldy*ldy;
1042
1043 Double_t cS = bq1*p2 + bq*bq1*(dy* px - dx* py) - bq*b;
1044 Double_t cS1= bq*p21 - bq*bq1*(dy*px1 - dx*py1) - bq1*b;
1045
1046 Double_t ca = bq*bq*bq1*d2 +2*( cS + bq*bq*(py1*dx-px1*dy)) ;
1047 Double_t ca1 = bq*bq1*bq1*d2 +2*( cS1 - bq1*bq1*(py*dx-px*dy)) ;
1048
1049 Double_t sa = 4*l2*p2 - ca*ca;
1050 Double_t sa1 = 4*l2*p21 - ca1*ca1;
e7b09c95 1051
f826d409 1052 if(sa<0) sa=0;
1053 if(sa1<0)sa1=0;
1054
1055 if( TMath::Abs(bq)>1.e-8){
1056 s = TMath::ATan2( bq*( bq1*(dx*px +dy*py) + a ) , cS )/bq;
1057 ds = TMath::ATan2(TMath::Sqrt(sa),ca)/bq;
1058 } else {
1059 s = ( (dx*px + dy*py) + (py*px1-px*py1)/bq1)/p2;
1060 ds = s*s - (d2-2*(px1*dy-py1*dx)/bq1)/p2;
1061 if( ds<0 ) ds = 0;
1062 ds = TMath::Sqrt(ds);
1063 }
1064
1065 if( TMath::Abs(bq1)>1.e-8){
1066 s1 = TMath::ATan2( -bq1*( bq*(dx*px1+dy*py1) + a), cS1 )/bq1;
1067 ds1 = TMath::ATan2(TMath::Sqrt(sa1),ca1)/bq1;
1068 } else {
1069 s1 = (-(dx*px1 + dy*py1) + (py*px1-px*py1)/bq)/p21;
1070 ds1 = s1*s1 - (d2+2*(px*dy-py*dx)/bq)/p21;
1071 if( ds1<0 ) ds1 = 0;
1072 ds1 = TMath::Sqrt(ds1);
1073 }
1074 }
1075
f826d409 1076 Double_t ss[2], ss1[2], g[2][5],g1[2][5];
1077
1078 ss[0] = s + ds;
1079 ss[1] = s - ds;
1080 ss1[0] = s1 + ds1;
1081 ss1[1] = s1 - ds1;
1082 for( Int_t i=0; i<2; i++){
1083 Double_t bs = bq*ss[i];
5fc72f28 1084 Double_t c = TMath::Cos(bs), sss = TMath::Sin(bs);
f826d409 1085 Double_t cB,sB;
1086 if( TMath::Abs(bq)>1.e-8){
1087 cB= (1-c)/bq;
5fc72f28 1088 sB= sss/bq;
f826d409 1089 }else{
bfd20868 1090 const Double_t kOvSqr6 = 1./TMath::Sqrt(6.);
1091 sB = (1.-bs*kOvSqr6)*(1.+bs*kOvSqr6)*ss[i];
f826d409 1092 cB = .5*sB*bs;
1093 }
1094 g[i][0] = fP[0] + sB*px + cB*py;
1095 g[i][1] = fP[1] - cB*px + sB*py;
1096 g[i][2] = fP[2] + ss[i]*pz;
5fc72f28 1097 g[i][3] = + c*px + sss*py;
1098 g[i][4] = - sss*px + c*py;
f826d409 1099
1100 bs = bq1*ss1[i];
5fc72f28 1101 c = TMath::Cos(bs); sss = TMath::Sin(bs);
f826d409 1102 if( TMath::Abs(bq1)>1.e-8){
1103 cB= (1-c)/bq1;
5fc72f28 1104 sB= sss/bq1;
f826d409 1105 }else{
bfd20868 1106 const Double_t kOvSqr6 = 1./TMath::Sqrt(6.);
1107 sB = (1.-bs*kOvSqr6)*(1.+bs*kOvSqr6)*ss1[i];
f826d409 1108 cB = .5*sB*bs;
1109 }
1110
e7b09c95 1111 g1[i][0] = p.fP[0] + sB*px1 + cB*py1;
1112 g1[i][1] = p.fP[1] - cB*px1 + sB*py1;
1113 g1[i][2] = p.fP[2] + ss[i]*pz1;
5fc72f28 1114 g1[i][3] = + c*px1 + sss*py1;
1115 g1[i][4] = - sss*px1 + c*py1;
f826d409 1116 }
1117
1118 Int_t i=0, i1=0;
1119
1120 Double_t dMin = 1.e10;
1121 for( Int_t j=0; j<2; j++){
1122 for( Int_t j1=0; j1<2; j1++){
1123 Double_t xx = g[j][0]-g1[j1][0];
1124 Double_t yy = g[j][1]-g1[j1][1];
1125 Double_t zz = g[j][2]-g1[j1][2];
1126 Double_t d = xx*xx + yy*yy + zz*zz;
1127 if( d<dMin ){
1128 dMin = d;
1129 i = j;
1130 i1 = j1;
1131 }
1132 }
1133 }
1134
1135 DS = ss[i];
1136 DS1 = ss1[i1];
91468b27 1137 if(0){
1138 Double_t x= g[i][0], y= g[i][1], z= g[i][2], ppx= g[i][3], ppy= g[i][4];
1139 Double_t x1=g1[i1][0], y1= g1[i1][1], z1= g1[i1][2], ppx1= g1[i1][3], ppy1= g1[i1][4];
1140 Double_t dx = x1-x;
1141 Double_t dy = y1-y;
1142 Double_t dz = z1-z;
1143 Double_t a = ppx*ppx1 + ppy*ppy1 + pz*pz1;
1144 Double_t b = dx*ppx1 + dy*ppy1 + dz*pz1;
1145 Double_t c = dx*ppx + dy*ppy + dz*pz ;
1146 Double_t pp2 = ppx*ppx + ppy*ppy + pz*pz;
1147 Double_t pp21= ppx1*ppx1 + ppy1*ppy1 + pz1*pz1;
706952f5 1148 Double_t det = pp2*pp21 - a*a;
91468b27 1149 if( TMath::Abs(det)>1.e-8 ){
1150 DS+=(a*b-pp21*c)/det;
1151 DS1+=(a*c-pp2*b)/det;
1152 }
f826d409 1153 }
1154}
1155
1156
1157
1158void AliKFParticleBase::TransportCBM( Double_t dS,
1159 Double_t P[], Double_t C[] ) const
1160{
1161 //* Transport the particle on dS, output to P[],C[], for CBM field
1162
1163 if( fQ==0 ){
1164 TransportLine( dS, P, C );
1165 return;
1166 }
1167
1168 const Double_t kCLight = 0.000299792458;
1169
1170 Double_t c = fQ*kCLight;
1171
1172 // construct coefficients
1173
1174 Double_t
1175 px = fP[3],
1176 py = fP[4],
1177 pz = fP[5];
1178
1179 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;
1180
1181 { // get field integrals
1182
1183 Double_t fld[3][3];
1184 Double_t p0[3], p1[3], p2[3];
1185
1186 // line track approximation
1187
1188 p0[0] = fP[0];
1189 p0[1] = fP[1];
1190 p0[2] = fP[2];
1191
1192 p2[0] = fP[0] + px*dS;
1193 p2[1] = fP[1] + py*dS;
1194 p2[2] = fP[2] + pz*dS;
1195
1196 p1[0] = 0.5*(p0[0]+p2[0]);
1197 p1[1] = 0.5*(p0[1]+p2[1]);
1198 p1[2] = 0.5*(p0[2]+p2[2]);
1199
1200 // first order track approximation
1201 {
1202 GetFieldValue( p0, fld[0] );
1203 GetFieldValue( p1, fld[1] );
1204 GetFieldValue( p2, fld[2] );
1205
1206 Double_t ssy1 = ( 7*fld[0][1] + 6*fld[1][1]-fld[2][1] )*c*dS*dS/96.;
1207 Double_t ssy2 = ( fld[0][1] + 2*fld[1][1] )*c*dS*dS/6.;
1208
1209 p1[0] -= ssy1*pz;
1210 p1[2] += ssy1*px;
1211 p2[0] -= ssy2*pz;
1212 p2[2] += ssy2*px;
1213 }
1214
1215 GetFieldValue( p0, fld[0] );
1216 GetFieldValue( p1, fld[1] );
1217 GetFieldValue( p2, fld[2] );
1218
1219 sx = c*( fld[0][0] + 4*fld[1][0] + fld[2][0] )*dS/6.;
1220 sy = c*( fld[0][1] + 4*fld[1][1] + fld[2][1] )*dS/6.;
1221 sz = c*( fld[0][2] + 4*fld[1][2] + fld[2][2] )*dS/6.;
1222
1223 ssx = c*( fld[0][0] + 2*fld[1][0])*dS*dS/6.;
1224 ssy = c*( fld[0][1] + 2*fld[1][1])*dS*dS/6.;
1225 ssz = c*( fld[0][2] + 2*fld[1][2])*dS*dS/6.;
1226
1227 Double_t c2[3][3] = { { 5, -4, -1},{ 44, 80, -4},{ 11, 44, 5} }; // /=360.
1228 Double_t cc2[3][3] = { { 38, 8, -4},{ 148, 208, -20},{ 3, 36, 3} }; // /=2520.
1229 for(Int_t n=0; n<3; n++)
1230 for(Int_t m=0; m<3; m++)
1231 {
1232 syz += c2[n][m]*fld[n][1]*fld[m][2];
1233 ssyz += cc2[n][m]*fld[n][1]*fld[m][2];
1234 }
1235
1236 syz *= c*c*dS*dS/360.;
1237 ssyz *= c*c*dS*dS*dS/2520.;
1238
1239 syy = c*( fld[0][1] + 4*fld[1][1] + fld[2][1] )*dS;
1240 syyy = syy*syy*syy / 1296;
1241 syy = syy*syy/72;
1242
1243 ssyy = ( fld[0][1]*( 38*fld[0][1] + 156*fld[1][1] - fld[2][1] )+
1244 fld[1][1]*( 208*fld[1][1] +16*fld[2][1] )+
1245 fld[2][1]*( 3*fld[2][1] )
1246 )*dS*dS*dS*c*c/2520.;
1247 ssyyy =
1248 (
1249 fld[0][1]*( fld[0][1]*( 85*fld[0][1] + 526*fld[1][1] - 7*fld[2][1] )+
1250 fld[1][1]*( 1376*fld[1][1] +84*fld[2][1] )+
1251 fld[2][1]*( 19*fld[2][1] ) )+
1252 fld[1][1]*( fld[1][1]*( 1376*fld[1][1] +256*fld[2][1] )+
1253 fld[2][1]*( 62*fld[2][1] ) )+
1254 fld[2][1]*fld[2][1] *( 3*fld[2][1] )
1255 )*dS*dS*dS*dS*c*c*c/90720.;
1256
1257 }
1258
1259 Double_t mJ[8][8];
1260 for( Int_t i=0; i<8; i++ ) for( Int_t j=0; j<8; j++) mJ[i][j]=0;
1261
1262 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;
1263 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;
1264 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;
1265
1266 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;
1267 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;
1268 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;
1269 mJ[6][6] = mJ[7][7] = 1;
1270
1271 P[0] = fP[0] + mJ[0][3]*px + mJ[0][4]*py + mJ[0][5]*pz;
1272 P[1] = fP[1] + mJ[1][3]*px + mJ[1][4]*py + mJ[1][5]*pz;
1273 P[2] = fP[2] + mJ[2][3]*px + mJ[2][4]*py + mJ[2][5]*pz;
1274 P[3] = mJ[3][3]*px + mJ[3][4]*py + mJ[3][5]*pz;
1275 P[4] = mJ[4][3]*px + mJ[4][4]*py + mJ[4][5]*pz;
1276 P[5] = mJ[5][3]*px + mJ[5][4]*py + mJ[5][5]*pz;
1277 P[6] = fP[6];
1278 P[7] = fP[7];
1279
1280 MultQSQt( mJ[0], fC, C);
1281
1282}
1283
1284
0d14829a 1285void AliKFParticleBase::TransportBz( Double_t b, Double_t t,
1286 Double_t p[], Double_t e[] ) const
f826d409 1287{
1288 //* Transport the particle on dS, output to P[],C[], for Bz field
616ffc76 1289
f826d409 1290 const Double_t kCLight = 0.000299792458;
5fc72f28 1291 b = b*fQ*kCLight;
0d14829a 1292 Double_t bs= b*t;
f826d409 1293 Double_t s = TMath::Sin(bs), c = TMath::Cos(bs);
1294 Double_t sB, cB;
616ffc76 1295 if( TMath::Abs(bs)>1.e-10){
5fc72f28 1296 sB= s/b;
1297 cB= (1-c)/b;
f826d409 1298 }else{
bfd20868 1299 const Double_t kOvSqr6 = 1./TMath::Sqrt(6.);
0d14829a 1300 sB = (1.-bs*kOvSqr6)*(1.+bs*kOvSqr6)*t;
f826d409 1301 cB = .5*sB*bs;
1302 }
1303
1304 Double_t px = fP[3];
1305 Double_t py = fP[4];
1306 Double_t pz = fP[5];
1307
5fc72f28 1308 p[0] = fP[0] + sB*px + cB*py;
1309 p[1] = fP[1] - cB*px + sB*py;
0d14829a 1310 p[2] = fP[2] + t*pz;
5fc72f28 1311 p[3] = c*px + s*py;
1312 p[4] = -s*px + c*py;
1313 p[5] = fP[5];
1314 p[6] = fP[6];
1315 p[7] = fP[7];
0d14829a 1316
1317 /*
f826d409 1318 Double_t mJ[8][8] = { {1,0,0, sB, cB, 0, 0, 0 },
616ffc76 1319 {0,1,0, -cB, sB, 0, 0, 0 },
0d14829a 1320 {0,0,1, 0, 0, t, 0, 0 },
616ffc76 1321 {0,0,0, c, s, 0, 0, 0 },
1322 {0,0,0, -s, c, 0, 0, 0 },
1323 {0,0,0, 0, 0, 1, 0, 0 },
1324 {0,0,0, 0, 0, 0, 1, 0 },
1325 {0,0,0, 0, 0, 0, 0, 1 } };
f826d409 1326 Double_t mA[8][8];
1327 for( Int_t k=0,i=0; i<8; i++)
1328 for( Int_t j=0; j<=i; j++, k++ ) mA[i][j] = mA[j][i] = fC[k];
1329
1330 Double_t mJC[8][8];
1331 for( Int_t i=0; i<8; i++ )
1332 for( Int_t j=0; j<8; j++ ){
1333 mJC[i][j]=0;
1334 for( Int_t k=0; k<8; k++ ) mJC[i][j]+=mJ[i][k]*mA[k][j];
1335 }
1336
1337 for( Int_t k=0,i=0; i<8; i++)
1338 for( Int_t j=0; j<=i; j++, k++ ){
0d14829a 1339 e[k] = 0;
1340 for( Int_t l=0; l<8; l++ ) e[k]+=mJC[i][l]*mJ[j][l];
f826d409 1341 }
1342
1343 return;
f826d409 1344 */
0d14829a 1345
1346 Double_t
1347 c6=fC[6], c7=fC[7], c8=fC[8], c17=fC[17], c18=fC[18],
1348 c24 = fC[24], c31 = fC[31];
1349
1350 Double_t
1351 cBC13 = cB*fC[13],
1352 mJC13 = c7 - cB*fC[9] + sB*fC[13],
1353 mJC14 = fC[11] - cBC13 + sB*fC[14],
1354 mJC23 = c8 + t*c18,
1355 mJC24 = fC[12] + t*fC[19],
1356 mJC33 = c*fC[9] + s*fC[13],
1357 mJC34 = c*fC[13] + s*fC[14],
1358 mJC43 = -s*fC[9] + c*fC[13],
1359 mJC44 = -s*fC[13] + c*fC[14];
1360
1361
1362 e[0]= fC[0] + 2*(sB*c6 + cB*fC[10]) + (sB*fC[9] + 2*cBC13)*sB + cB*cB*fC[14];
1363 e[1]= fC[1] - cB*c6 + sB*fC[10] + mJC13*sB + mJC14*cB;
1364 e[2]= fC[2] - cB*c7 + sB*fC[11] - mJC13*cB + mJC14*sB;
1365 e[3]= fC[3] + t*fC[15] + mJC23*sB + mJC24*cB;
1366 e[4]= fC[4] + t*fC[16] - mJC23*cB + mJC24*sB;
1367
1368 e[15]= fC[15] + c18*sB + fC[19]*cB;
1369 e[16]= fC[16] - c18*cB + fC[19]*sB;
1370 e[17]= c17 + fC[20]*t;
1371 e[18]= c18*c + fC[19]*s;
1372 e[19]= -c18*s + fC[19]*c;
1373
1374 e[5]= fC[5] + (c17 + e[17] )*t;
1375
1376 e[6]= c*c6 + s*fC[10] + mJC33*sB + mJC34*cB;
1377 e[7]= c*c7 + s*fC[11] - mJC33*cB + mJC34*sB;
1378 e[8]= c*c8 + s*fC[12] + e[18]*t;
1379 e[9]= mJC33*c + mJC34*s;
1380 e[10]= -s*c6 + c*fC[10] + mJC43*sB + mJC44*cB;
1381
1382
1383 e[11]= -s*c7 + c*fC[11] - mJC43*cB + mJC44*sB;
1384 e[12]= -s*c8 + c*fC[12] + e[19]*t;
1385 e[13]= mJC43*c + mJC44*s;
1386 e[14]= -mJC43*s + mJC44*c;
1387 e[20]= fC[20];
1388 e[21]= fC[21] + fC[25]*cB + c24*sB;
1389 e[22]= fC[22] - c24*cB + fC[25]*sB;
1390 e[23]= fC[23] + fC[26]*t;
1391 e[24]= c*c24 + s*fC[25];
1392 e[25]= c*fC[25] - c24*s;
1393 e[26]= fC[26];
1394 e[27]= fC[27];
1395 e[28]= fC[28] + fC[32]*cB + c31*sB;
1396 e[29]= fC[29] - c31*cB + fC[32]*sB;
1397 e[30]= fC[30] + fC[33]*t;
1398 e[31]= c*c31 + s*fC[32];
1399 e[32]= c*fC[32] - s*c31;
1400 e[33]= fC[33];
1401 e[34]= fC[34];
1402 e[35]= fC[35];
f826d409 1403}
1404
1405
1406Double_t AliKFParticleBase::GetDistanceFromVertex( const AliKFParticleBase &Vtx ) const
1407{
1408 //* Calculate distance from vertex [cm]
1409
1410 return GetDistanceFromVertex( Vtx.fP );
1411}
1412
1413Double_t AliKFParticleBase::GetDistanceFromVertex( const Double_t vtx[] ) const
1414{
1415 //* Calculate distance from vertex [cm]
1416
1417 Double_t mP[8], mC[36];
1418 Transport( GetDStoPoint(vtx), mP, mC );
1419 Double_t d[3]={ vtx[0]-mP[0], vtx[1]-mP[1], vtx[2]-mP[2]};
1420 return TMath::Sqrt( d[0]*d[0]+d[1]*d[1]+d[2]*d[2] );
1421}
1422
616ffc76 1423Double_t AliKFParticleBase::GetDistanceFromParticle( const AliKFParticleBase &p )
f826d409 1424 const
1425{
616ffc76 1426 //* Calculate distance to other particle [cm]
f826d409 1427
1428 Double_t dS, dS1;
616ffc76 1429 GetDStoParticle( p, dS, dS1 );
f826d409 1430 Double_t mP[8], mC[36], mP1[8], mC1[36];
616ffc76 1431 Transport( dS, mP, mC );
1432 p.Transport( dS1, mP1, mC1 );
f826d409 1433 Double_t dx = mP[0]-mP1[0];
1434 Double_t dy = mP[1]-mP1[1];
1435 Double_t dz = mP[2]-mP1[2];
55ac3e1e 1436 dz = 0;
f826d409 1437 return TMath::Sqrt(dx*dx+dy*dy+dz*dz);
1438}
1439
1440Double_t AliKFParticleBase::GetDeviationFromVertex( const AliKFParticleBase &Vtx ) const
1441{
1442 //* Calculate sqrt(Chi2/ndf) deviation from vertex
1443
1444 return GetDeviationFromVertex( Vtx.fP, Vtx.fC );
1445}
1446
1447
1448Double_t AliKFParticleBase::GetDeviationFromVertex( const Double_t v[], const Double_t Cv[] ) const
1449{
1450 //* Calculate sqrt(Chi2/ndf) deviation from vertex
1451 //* v = [xyz], Cv=[Cxx,Cxy,Cyy,Cxz,Cyz,Czz]-covariance matrix
1452
1453 Double_t mP[8];
1454 Double_t mC[36];
1455
1456 Transport( GetDStoPoint(v), mP, mC );
1457
1458 Double_t d[3]={ v[0]-mP[0], v[1]-mP[1], v[2]-mP[2]};
1459
1460 Double_t sigmaS = .1+10.*TMath::Sqrt( (d[0]*d[0]+d[1]*d[1]+d[2]*d[2])/
1461 (mP[3]*mP[3]+mP[4]*mP[4]+mP[5]*mP[5]) );
1462
1463
1464 Double_t h[3] = { mP[3]*sigmaS, mP[4]*sigmaS, mP[5]*sigmaS };
1465
1466 Double_t mSi[6] =
1467 { mC[0] +h[0]*h[0],
1468 mC[1] +h[1]*h[0], mC[2] +h[1]*h[1],
1469 mC[3] +h[2]*h[0], mC[4] +h[2]*h[1], mC[5] +h[2]*h[2] };
1470
1471 if( Cv ){
1472 mSi[0]+=Cv[0];
1473 mSi[1]+=Cv[1];
1474 mSi[2]+=Cv[2];
1475 mSi[3]+=Cv[3];
1476 mSi[4]+=Cv[4];
1477 mSi[5]+=Cv[5];
1478 }
1479
1480 Double_t mS[6];
1481
1482 mS[0] = mSi[2]*mSi[5] - mSi[4]*mSi[4];
1483 mS[1] = mSi[3]*mSi[4] - mSi[1]*mSi[5];
1484 mS[2] = mSi[0]*mSi[5] - mSi[3]*mSi[3];
1485 mS[3] = mSi[1]*mSi[4] - mSi[2]*mSi[3];
1486 mS[4] = mSi[1]*mSi[3] - mSi[0]*mSi[4];
1487 mS[5] = mSi[0]*mSi[2] - mSi[1]*mSi[1];
1488
1489 Double_t s = ( mSi[0]*mS[0] + mSi[1]*mS[1] + mSi[3]*mS[3] );
1490 s = ( s > 1.E-20 ) ?1./s :0;
1491
1492 return TMath::Sqrt( TMath::Abs(s*( ( mS[0]*d[0] + mS[1]*d[1] + mS[3]*d[2])*d[0]
1493 +(mS[1]*d[0] + mS[2]*d[1] + mS[4]*d[2])*d[1]
1494 +(mS[3]*d[0] + mS[4]*d[1] + mS[5]*d[2])*d[2] ))/2);
1495}
1496
1497
616ffc76 1498Double_t AliKFParticleBase::GetDeviationFromParticle( const AliKFParticleBase &p )
f826d409 1499 const
1500{
616ffc76 1501 //* Calculate sqrt(Chi2/ndf) deviation from other particle
f826d409 1502
1503 Double_t dS, dS1;
616ffc76 1504 GetDStoParticle( p, dS, dS1 );
f826d409 1505 Double_t mP1[8], mC1[36];
616ffc76 1506 p.Transport( dS1, mP1, mC1 );
f826d409 1507
1508 Double_t d[3]={ fP[0]-mP1[0], fP[1]-mP1[1], fP[2]-mP1[2]};
1509
1510 Double_t sigmaS = .1+10.*TMath::Sqrt( (d[0]*d[0]+d[1]*d[1]+d[2]*d[2])/
616ffc76 1511 (mP1[3]*mP1[3]+mP1[4]*mP1[4]+mP1[5]*mP1[5]) );
f826d409 1512
1513 Double_t h[3] = { mP1[3]*sigmaS, mP1[4]*sigmaS, mP1[5]*sigmaS };
1514
1515 mC1[0] +=h[0]*h[0];
1516 mC1[1] +=h[1]*h[0];
1517 mC1[2] +=h[1]*h[1];
1518 mC1[3] +=h[2]*h[0];
1519 mC1[4] +=h[2]*h[1];
1520 mC1[5] +=h[2]*h[2];
1521
1522 return GetDeviationFromVertex( mP1, mC1 )*TMath::Sqrt(2./1.);
1523}
1524
1525
1526
de0d0ceb 1527void AliKFParticleBase::SubtractFromVertex( AliKFParticleBase &Vtx ) const
f826d409 1528{
1529 //* Subtract the particle from the vertex
de0d0ceb 1530
f826d409 1531 Double_t fld[3];
1532 {
de0d0ceb 1533 GetFieldValue( Vtx.fP, fld );
f826d409 1534 const Double_t kCLight = 0.000299792458;
1535 fld[0]*=kCLight; fld[1]*=kCLight; fld[2]*=kCLight;
1536 }
1537
1538 Double_t m[8];
1539 Double_t mCm[36];
f826d409 1540
de0d0ceb 1541 if( Vtx.fIsLinearized ){
1542 GetMeasurement( Vtx.fVtxGuess, m, mCm );
1543 } else {
1544 GetMeasurement( Vtx.fP, m, mCm );
f826d409 1545 }
de0d0ceb 1546
f826d409 1547 Double_t mV[6];
1548
de0d0ceb 1549 mV[ 0] = mCm[ 0];
1550 mV[ 1] = mCm[ 1];
1551 mV[ 2] = mCm[ 2];
1552 mV[ 3] = mCm[ 3];
1553 mV[ 4] = mCm[ 4];
1554 mV[ 5] = mCm[ 5];
f826d409 1555
1556 //*
1557
1558 Double_t mS[6];
1559 {
de0d0ceb 1560 Double_t mSi[6] = { mV[0]-Vtx.fC[0],
1561 mV[1]-Vtx.fC[1], mV[2]-Vtx.fC[2],
1562 mV[3]-Vtx.fC[3], mV[4]-Vtx.fC[4], mV[5]-Vtx.fC[5] };
f826d409 1563
1564 mS[0] = mSi[2]*mSi[5] - mSi[4]*mSi[4];
1565 mS[1] = mSi[3]*mSi[4] - mSi[1]*mSi[5];
1566 mS[2] = mSi[0]*mSi[5] - mSi[3]*mSi[3];
1567 mS[3] = mSi[1]*mSi[4] - mSi[2]*mSi[3];
1568 mS[4] = mSi[1]*mSi[3] - mSi[0]*mSi[4];
1569 mS[5] = mSi[0]*mSi[2] - mSi[1]*mSi[1];
1570
1571 Double_t s = ( mSi[0]*mS[0] + mSi[1]*mS[1] + mSi[3]*mS[3] );
1572 s = ( s > 1.E-20 ) ?1./s :0;
1573 mS[0]*=s;
1574 mS[1]*=s;
1575 mS[2]*=s;
1576 mS[3]*=s;
1577 mS[4]*=s;
1578 mS[5]*=s;
1579 }
1580
1581 //* Residual (measured - estimated)
1582
de0d0ceb 1583 Double_t zeta[3] = { m[0]-Vtx.fP[0], m[1]-Vtx.fP[1], m[2]-Vtx.fP[2] };
f826d409 1584
1585 //* mCHt = mCH' - D'
1586
1587 Double_t mCHt0[3], mCHt1[3], mCHt2[3];
1588
de0d0ceb 1589 mCHt0[0]=Vtx.fC[ 0] ; mCHt1[0]=Vtx.fC[ 1] ; mCHt2[0]=Vtx.fC[ 3] ;
1590 mCHt0[1]=Vtx.fC[ 1] ; mCHt1[1]=Vtx.fC[ 2] ; mCHt2[1]=Vtx.fC[ 4] ;
1591 mCHt0[2]=Vtx.fC[ 3] ; mCHt1[2]=Vtx.fC[ 4] ; mCHt2[2]=Vtx.fC[ 5] ;
f826d409 1592
1593 //* Kalman gain K = mCH'*S
1594
1595 Double_t k0[3], k1[3], k2[3];
1596
1597 for(Int_t i=0;i<3;++i){
1598 k0[i] = mCHt0[i]*mS[0] + mCHt1[i]*mS[1] + mCHt2[i]*mS[3];
1599 k1[i] = mCHt0[i]*mS[1] + mCHt1[i]*mS[2] + mCHt2[i]*mS[4];
1600 k2[i] = mCHt0[i]*mS[3] + mCHt1[i]*mS[4] + mCHt2[i]*mS[5];
1601 }
1602
1603 //* New estimation of the vertex position r += K*zeta
1604
de0d0ceb 1605 Double_t dChi2 = -(mS[0]*zeta[0] + mS[1]*zeta[1] + mS[3]*zeta[2])*zeta[0]
1606 + (mS[1]*zeta[0] + mS[2]*zeta[1] + mS[4]*zeta[2])*zeta[1]
1607 + (mS[3]*zeta[0] + mS[4]*zeta[1] + mS[5]*zeta[2])*zeta[2];
1608
1609 if( Vtx.fChi2 - dChi2 < 0 ) return;
1610
f826d409 1611 for(Int_t i=0;i<3;++i)
de0d0ceb 1612 Vtx.fP[i] -= k0[i]*zeta[0] + k1[i]*zeta[1] + k2[i]*zeta[2];
f826d409 1613
1614 //* New covariance matrix C -= K*(mCH')'
1615
1616 for(Int_t i=0, k=0;i<3;++i){
1617 for(Int_t j=0;j<=i;++j,++k)
de0d0ceb 1618 Vtx.fC[k] += k0[i]*mCHt0[j] + k1[i]*mCHt1[j] + k2[i]*mCHt2[j];
f826d409 1619 }
1620
1621 //* Calculate Chi^2
de0d0ceb 1622
1623 Vtx.fNDF -= 2;
1624 Vtx.fChi2 -= dChi2;
f826d409 1625}
1626
1627
1628
1629void AliKFParticleBase::TransportLine( Double_t dS,
1630 Double_t P[], Double_t C[] ) const
1631{
1632 //* Transport the particle as a straight line
1633
1634 P[0] = fP[0] + dS*fP[3];
1635 P[1] = fP[1] + dS*fP[4];
1636 P[2] = fP[2] + dS*fP[5];
1637 P[3] = fP[3];
1638 P[4] = fP[4];
1639 P[5] = fP[5];
1640 P[6] = fP[6];
1641 P[7] = fP[7];
1642
1643 Double_t c6 = fC[ 6] + dS*fC[ 9];
1644 Double_t c11 = fC[11] + dS*fC[14];
1645 Double_t c17 = fC[17] + dS*fC[20];
1646 Double_t sc13 = dS*fC[13];
1647 Double_t sc18 = dS*fC[18];
1648 Double_t sc19 = dS*fC[19];
1649
1650 C[ 0] = fC[ 0] + dS*( fC[ 6] + c6 );
1651 C[ 2] = fC[ 2] + dS*( fC[11] + c11 );
1652 C[ 5] = fC[ 5] + dS*( fC[17] + c17 );
1653
1654 C[ 7] = fC[ 7] + sc13;
1655 C[ 8] = fC[ 8] + sc18;
1656 C[ 9] = fC[ 9];
1657
1658 C[12] = fC[12] + sc19;
1659
1660 C[ 1] = fC[ 1] + dS*( fC[10] + C[ 7] );
1661 C[ 3] = fC[ 3] + dS*( fC[15] + C[ 8] );
1662 C[ 4] = fC[ 4] + dS*( fC[16] + C[12] );
1663 C[ 6] = c6;
1664
1665 C[10] = fC[10] + sc13;
1666 C[11] = c11;
1667
1668 C[13] = fC[13];
1669 C[14] = fC[14];
1670 C[15] = fC[15] + sc18;
1671 C[16] = fC[16] + sc19;
1672 C[17] = c17;
1673
1674 C[18] = fC[18];
1675 C[19] = fC[19];
1676 C[20] = fC[20];
1677 C[21] = fC[21] + dS*fC[24];
1678 C[22] = fC[22] + dS*fC[25];
1679 C[23] = fC[23] + dS*fC[26];
1680
1681 C[24] = fC[24];
1682 C[25] = fC[25];
1683 C[26] = fC[26];
1684 C[27] = fC[27];
1685 C[28] = fC[28] + dS*fC[31];
1686 C[29] = fC[29] + dS*fC[32];
1687 C[30] = fC[30] + dS*fC[33];
1688
1689 C[31] = fC[31];
1690 C[32] = fC[32];
1691 C[33] = fC[33];
1692 C[34] = fC[34];
1693 C[35] = fC[35];
1694}
1695
1696
a65041d0 1697void AliKFParticleBase::ConstructGammaBz( const AliKFParticleBase &daughter1,
1698 const AliKFParticleBase &daughter2, double Bz )
1699{
1700 //* Create gamma
1701
1702 const AliKFParticleBase *daughters[2] = { &daughter1, &daughter2};
1703
1704 double v0[3];
1705
1706 if( !fIsLinearized ){
1707 Double_t ds, ds1;
1708 Double_t m[8];
1709 Double_t mCd[36];
1710 daughter1.GetDStoParticle(daughter2, ds, ds1);
1711 daughter1.Transport( ds, m, mCd );
55ac3e1e 1712 fP[0] = m[0];
1713 fP[1] = m[1];
1714 fP[2] = m[2];
a65041d0 1715 daughter2.Transport( ds1, m, mCd );
55ac3e1e 1716 fP[0] = .5*( fP[0] + m[0] );
1717 fP[1] = .5*( fP[1] + m[1] );
1718 fP[2] = .5*( fP[2] + m[2] );
a65041d0 1719 } else {
55ac3e1e 1720 fP[0] = fVtxGuess[0];
1721 fP[1] = fVtxGuess[1];
1722 fP[2] = fVtxGuess[2];
a65041d0 1723 }
1724
a65041d0 1725 double daughterP[2][8], daughterC[2][36];
1726 double vtxMom[2][3];
1727
55ac3e1e 1728 int nIter = fIsLinearized ?1 :2;
a65041d0 1729
55ac3e1e 1730 for( int iter=0; iter<nIter; iter++){
a65041d0 1731
55ac3e1e 1732 v0[0] = fP[0];
1733 v0[1] = fP[1];
1734 v0[2] = fP[2];
a65041d0 1735
55ac3e1e 1736 fAtProductionVertex = 0;
1737 fSFromDecay = 0;
1738 fP[0] = v0[0];
1739 fP[1] = v0[1];
1740 fP[2] = v0[2];
1741 fP[3] = 0;
1742 fP[4] = 0;
1743 fP[5] = 0;
1744 fP[6] = 0;
1745 fP[7] = 0;
1746
1747
1748 // fit daughters to the vertex guess
a65041d0 1749
55ac3e1e 1750 {
1751 for( int id=0; id<2; id++ ){
1752
1753 double *p = daughterP[id];
1754 double *mC = daughterC[id];
1755
1756 daughters[id]->GetMeasurement( v0, p, mC );
1757
1758 Double_t mAi[6];
1759 InvertSym3(mC, mAi );
1760
1761 Double_t mB[3][3];
1762
1763 mB[0][0] = mC[ 6]*mAi[0] + mC[ 7]*mAi[1] + mC[ 8]*mAi[3];
1764 mB[0][1] = mC[ 6]*mAi[1] + mC[ 7]*mAi[2] + mC[ 8]*mAi[4];
1765 mB[0][2] = mC[ 6]*mAi[3] + mC[ 7]*mAi[4] + mC[ 8]*mAi[5];
1766
1767 mB[1][0] = mC[10]*mAi[0] + mC[11]*mAi[1] + mC[12]*mAi[3];
1768 mB[1][1] = mC[10]*mAi[1] + mC[11]*mAi[2] + mC[12]*mAi[4];
1769 mB[1][2] = mC[10]*mAi[3] + mC[11]*mAi[4] + mC[12]*mAi[5];
1770
1771 mB[2][0] = mC[15]*mAi[0] + mC[16]*mAi[1] + mC[17]*mAi[3];
1772 mB[2][1] = mC[15]*mAi[1] + mC[16]*mAi[2] + mC[17]*mAi[4];
1773 mB[2][2] = mC[15]*mAi[3] + mC[16]*mAi[4] + mC[17]*mAi[5];
1774
1775 Double_t z[3] = { v0[0]-p[0], v0[1]-p[1], v0[2]-p[2] };
1776
1777 vtxMom[id][0] = p[3] + mB[0][0]*z[0] + mB[0][1]*z[1] + mB[0][2]*z[2];
1778 vtxMom[id][1] = p[4] + mB[1][0]*z[0] + mB[1][1]*z[1] + mB[1][2]*z[2];
1779 vtxMom[id][2] = p[5] + mB[2][0]*z[0] + mB[2][1]*z[1] + mB[2][2]*z[2];
1780
1781 daughters[id]->Transport( daughters[id]->GetDStoPoint(v0), p, mC );
a65041d0 1782
55ac3e1e 1783 }
a65041d0 1784
55ac3e1e 1785 } // fit daughters to guess
1786
1787
1788 // fit new vertex
1789 {
a65041d0 1790
55ac3e1e 1791 double mpx0 = vtxMom[0][0]+vtxMom[1][0];
1792 double mpy0 = vtxMom[0][1]+vtxMom[1][1];
1793 double mpt0 = TMath::Sqrt(mpx0*mpx0 + mpy0*mpy0);
1794 // double a0 = TMath::ATan2(mpy0,mpx0);
a65041d0 1795
55ac3e1e 1796 double ca0 = mpx0/mpt0;
1797 double sa0 = mpy0/mpt0;
1798 double r[3] = { v0[0], v0[1], v0[2] };
1799 double mC[3][3] = {{10000., 0 , 0 },
1800 {0, 10000., 0 },
1801 {0, 0, 10000. } };
1802 double chi2=0;
a65041d0 1803
55ac3e1e 1804 for( int id=0; id<2; id++ ){
1805 const Double_t kCLight = 0.000299792458;
1806 Double_t q = Bz*daughters[id]->GetQ()*kCLight;
1807 Double_t px0 = vtxMom[id][0];
1808 Double_t py0 = vtxMom[id][1];
1809 Double_t pz0 = vtxMom[id][2];
1810 Double_t pt0 = TMath::Sqrt(px0*px0+py0*py0);
1811 Double_t mG[3][6], mB[3], mH[3][3];
1812 // r = {vx,vy,vz};
1813 // m = {x,y,z,Px,Py,Pz};
1814 // V = daughter.C
1815 // G*m + B = H*r;
1816 // q*x + Py - q*vx - sin(a)*Pt = 0
1817 // q*y - Px - q*vy + cos(a)*Pt = 0
1818 // (Px*cos(a) + Py*sin(a) ) (vz -z) - Pz( cos(a)*(vx-x) + sin(a)*(vy-y)) = 0
1819
1820 mG[0][0] = q;
1821 mG[0][1] = 0;
1822 mG[0][2] = 0;
1823 mG[0][3] = -sa0*px0/pt0;
1824 mG[0][4] = 1 -sa0*py0/pt0;
1825 mG[0][5] = 0;
1826 mH[0][0] = q;
1827 mH[0][1] = 0;
1828 mH[0][2] = 0;
1829 mB[0] = py0 - sa0*pt0 - mG[0][3]*px0 - mG[0][4]*py0 ;
1830
1831 // q*y - Px - q*vy + cos(a)*Pt = 0
1832
1833 mG[1][0] = 0;
1834 mG[1][1] = q;
1835 mG[1][2] = 0;
1836 mG[1][3] = -1 + ca0*px0/pt0;
1837 mG[1][4] = + ca0*py0/pt0;
1838 mG[1][5] = 0;
1839 mH[1][0] = 0;
1840 mH[1][1] = q;
1841 mH[1][2] = 0;
1842 mB[1] = -px0 + ca0*pt0 - mG[1][3]*px0 - mG[1][4]*py0 ;
1843
1844 // (Px*cos(a) + Py*sin(a) ) (z -vz) - Pz( cos(a)*(x-vx) + sin(a)*(y-vy)) = 0
a65041d0 1845
55ac3e1e 1846 mG[2][0] = -pz0*ca0;
1847 mG[2][1] = -pz0*sa0;
1848 mG[2][2] = px0*ca0 + py0*sa0;
1849 mG[2][3] = 0;
1850 mG[2][4] = 0;
1851 mG[2][5] = 0;
1852
1853 mH[2][0] = mG[2][0];
1854 mH[2][1] = mG[2][1];
1855 mH[2][2] = mG[2][2];
1856
1857 mB[2] = 0;
1858
1859 // fit the vertex
a65041d0 1860
55ac3e1e 1861 // V = GVGt
a65041d0 1862
55ac3e1e 1863 double mGV[3][6];
1864 double mV[6];
1865 double m[3];
1866 for( int i=0; i<3; i++ ){
1867 m[i] = mB[i];
1868 for( int k=0; k<6; k++ ) m[i]+=mG[i][k]*daughterP[id][k];
a65041d0 1869 }
55ac3e1e 1870 for( int i=0; i<3; i++ ){
1871 for( int j=0; j<6; j++ ){
1872 mGV[i][j] = 0;
1873 for( int k=0; k<6; k++ ) mGV[i][j]+=mG[i][k]*daughterC[id][ IJ(k,j) ];
1874 }
a65041d0 1875 }
55ac3e1e 1876 for( int i=0, k=0; i<3; i++ ){
1877 for( int j=0; j<=i; j++,k++ ){
1878 mV[k] = 0;
1879 for( int l=0; l<6; l++ ) mV[k]+=mGV[i][l]*mG[j][l];
1880 }
1881 }
1882
a65041d0 1883
55ac3e1e 1884 //* CHt
1885
1886 Double_t mCHt[3][3];
1887 Double_t mHCHt[6];
1888 Double_t mHr[3];
1889 for( int i=0; i<3; i++ ){
1890 mHr[i] = 0;
1891 for( int k=0; k<3; k++ ) mHr[i]+= mH[i][k]*r[k];
1892 }
a65041d0 1893
55ac3e1e 1894 for( int i=0; i<3; i++ ){
1895 for( int j=0; j<3; j++){
1896 mCHt[i][j] = 0;
1897 for( int k=0; k<3; k++ ) mCHt[i][j]+= mC[i][k]*mH[j][k];
1898 }
a65041d0 1899 }
a65041d0 1900
55ac3e1e 1901 for( int i=0, k=0; i<3; i++ ){
1902 for( int j=0; j<=i; j++, k++ ){
1903 mHCHt[k] = 0;
1904 for( int l=0; l<3; l++ ) mHCHt[k]+= mH[i][l]*mCHt[l][j];
1905 }
a65041d0 1906 }
a65041d0 1907
55ac3e1e 1908 Double_t mS[6] = { mHCHt[0]+mV[0],
1909 mHCHt[1]+mV[1], mHCHt[2]+mV[2],
1910 mHCHt[3]+mV[3], mHCHt[4]+mV[4], mHCHt[5]+mV[5] };
a65041d0 1911
1912
55ac3e1e 1913 InvertSym3(mS,mS);
a65041d0 1914
55ac3e1e 1915 //* Residual (measured - estimated)
a65041d0 1916
55ac3e1e 1917 Double_t zeta[3] = { m[0]-mHr[0], m[1]-mHr[1], m[2]-mHr[2] };
a65041d0 1918
55ac3e1e 1919 //* Kalman gain K = mCH'*S
a65041d0 1920
55ac3e1e 1921 Double_t k[3][3];
a65041d0 1922
55ac3e1e 1923 for(Int_t i=0;i<3;++i){
1924 k[i][0] = mCHt[i][0]*mS[0] + mCHt[i][1]*mS[1] + mCHt[i][2]*mS[3];
1925 k[i][1] = mCHt[i][0]*mS[1] + mCHt[i][1]*mS[2] + mCHt[i][2]*mS[4];
1926 k[i][2] = mCHt[i][0]*mS[3] + mCHt[i][1]*mS[4] + mCHt[i][2]*mS[5];
1927 }
a65041d0 1928
55ac3e1e 1929 //* New estimation of the vertex position r += K*zeta
a65041d0 1930
55ac3e1e 1931 for(Int_t i=0;i<3;++i)
1932 r[i] = r[i] + k[i][0]*zeta[0] + k[i][1]*zeta[1] + k[i][2]*zeta[2];
a65041d0 1933
55ac3e1e 1934 //* New covariance matrix C -= K*(mCH')'
a65041d0 1935
55ac3e1e 1936 for(Int_t i=0;i<3;++i){
1937 for(Int_t j=0;j<=i;++j){
1938 mC[i][j] = mC[i][j] - (k[i][0]*mCHt[j][0] + k[i][1]*mCHt[j][1] + k[i][2]*mCHt[j][2]);
1939 mC[j][i] = mC[i][j];
1940 }
a65041d0 1941 }
a65041d0 1942
1943
55ac3e1e 1944 //* Calculate Chi^2
a65041d0 1945
55ac3e1e 1946 chi2 += ( ( mS[0]*zeta[0] + mS[1]*zeta[1] + mS[3]*zeta[2] )*zeta[0]
1947 +(mS[1]*zeta[0] + mS[2]*zeta[1] + mS[4]*zeta[2] )*zeta[1]
1948 +(mS[3]*zeta[0] + mS[4]*zeta[1] + mS[5]*zeta[2] )*zeta[2] );
1949 }
a65041d0 1950
55ac3e1e 1951 // store vertex
a65041d0 1952
55ac3e1e 1953 fNDF = 2;
1954 fChi2 = chi2;
1955 for( int i=0; i<3; i++ ) fP[i] = r[i];
1956 for( int i=0,k=0; i<3; i++ ){
1957 for( int j=0; j<=i; j++,k++ ){
1958 fC[k] = mC[i][j];
1959 }
a65041d0 1960 }
1961 }
55ac3e1e 1962
1963 } // iterations
a65041d0 1964
1965 // now fit daughters to the vertex
1966
1967 fQ = 0;
1968 fSFromDecay = 0;
1969
1970 for(Int_t i=3;i<8;++i) fP[i]=0.;
55ac3e1e 1971 for(Int_t i=6;i<35;++i) fC[i]=0.;
1972 fC[35] = 100.;
a65041d0 1973
1974 for( int id=0; id<2; id++ ){
1975
1976 double *p = daughterP[id];
1977 double *mC = daughterC[id];
1978 daughters[id]->GetMeasurement( v0, p, mC );
1979
1980 const Double_t *m = fP, *mV = fC;
1981
1982 Double_t mAi[6];
1983 InvertSym3(mC, mAi );
1984
1985 Double_t mB[4][3];
1986
1987 mB[0][0] = mC[ 6]*mAi[0] + mC[ 7]*mAi[1] + mC[ 8]*mAi[3];
1988 mB[0][1] = mC[ 6]*mAi[1] + mC[ 7]*mAi[2] + mC[ 8]*mAi[4];
1989 mB[0][2] = mC[ 6]*mAi[3] + mC[ 7]*mAi[4] + mC[ 8]*mAi[5];
1990
1991 mB[1][0] = mC[10]*mAi[0] + mC[11]*mAi[1] + mC[12]*mAi[3];
1992 mB[1][1] = mC[10]*mAi[1] + mC[11]*mAi[2] + mC[12]*mAi[4];
1993 mB[1][2] = mC[10]*mAi[3] + mC[11]*mAi[4] + mC[12]*mAi[5];
1994
1995 mB[2][0] = mC[15]*mAi[0] + mC[16]*mAi[1] + mC[17]*mAi[3];
1996 mB[2][1] = mC[15]*mAi[1] + mC[16]*mAi[2] + mC[17]*mAi[4];
1997 mB[2][2] = mC[15]*mAi[3] + mC[16]*mAi[4] + mC[17]*mAi[5];
1998
1999 mB[3][0] = mC[21]*mAi[0] + mC[22]*mAi[1] + mC[23]*mAi[3];
2000 mB[3][1] = mC[21]*mAi[1] + mC[22]*mAi[2] + mC[23]*mAi[4];
2001 mB[3][2] = mC[21]*mAi[3] + mC[22]*mAi[4] + mC[23]*mAi[5];
2002
2003
2004 Double_t z[3] = { m[0]-p[0], m[1]-p[1], m[2]-p[2] };
2005
2006 {
2007 Double_t mAV[6] = { mC[0]-mV[0], mC[1]-mV[1], mC[2]-mV[2],
2008 mC[3]-mV[3], mC[4]-mV[4], mC[5]-mV[5] };
2009
2010 Double_t mAVi[6];
2011 if( !InvertSym3(mAV, mAVi) ){
2012 Double_t dChi2 = ( +(mAVi[0]*z[0] + mAVi[1]*z[1] + mAVi[3]*z[2])*z[0]
2013 +(mAVi[1]*z[0] + mAVi[2]*z[1] + mAVi[4]*z[2])*z[1]
2014 +(mAVi[3]*z[0] + mAVi[4]*z[1] + mAVi[5]*z[2])*z[2] );
2015 fChi2+= TMath::Abs( dChi2 );
2016 }
2017 fNDF += 2;
2018 }
2019
2020 //* Add the daughter momentum to the particle momentum
2021
2022 fP[3]+= p[3] + mB[0][0]*z[0] + mB[0][1]*z[1] + mB[0][2]*z[2];
2023 fP[4]+= p[4] + mB[1][0]*z[0] + mB[1][1]*z[1] + mB[1][2]*z[2];
2024 fP[5]+= p[5] + mB[2][0]*z[0] + mB[2][1]*z[1] + mB[2][2]*z[2];
2025 fP[6]+= p[6] + mB[3][0]*z[0] + mB[3][1]*z[1] + mB[3][2]*z[2];
2026
2027 Double_t d0, d1, d2;
2028
2029 d0= mB[0][0]*mV[0] + mB[0][1]*mV[1] + mB[0][2]*mV[3] - mC[ 6];
2030 d1= mB[0][0]*mV[1] + mB[0][1]*mV[2] + mB[0][2]*mV[4] - mC[ 7];
2031 d2= mB[0][0]*mV[3] + mB[0][1]*mV[4] + mB[0][2]*mV[5] - mC[ 8];
2032
2033 //fC[6]+= mC[ 6] + d0;
2034 //fC[7]+= mC[ 7] + d1;
2035 //fC[8]+= mC[ 8] + d2;
2036 fC[9]+= mC[ 9] + d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
2037
2038 d0= mB[1][0]*mV[0] + mB[1][1]*mV[1] + mB[1][2]*mV[3] - mC[10];
2039 d1= mB[1][0]*mV[1] + mB[1][1]*mV[2] + mB[1][2]*mV[4] - mC[11];
2040 d2= mB[1][0]*mV[3] + mB[1][1]*mV[4] + mB[1][2]*mV[5] - mC[12];
2041
2042 //fC[10]+= mC[10]+ d0;
2043 //fC[11]+= mC[11]+ d1;
2044 //fC[12]+= mC[12]+ d2;
2045 fC[13]+= mC[13]+ d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
2046 fC[14]+= mC[14]+ d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
2047
2048 d0= mB[2][0]*mV[0] + mB[2][1]*mV[1] + mB[2][2]*mV[3] - mC[15];
2049 d1= mB[2][0]*mV[1] + mB[2][1]*mV[2] + mB[2][2]*mV[4] - mC[16];
2050 d2= mB[2][0]*mV[3] + mB[2][1]*mV[4] + mB[2][2]*mV[5] - mC[17];
2051
2052 //fC[15]+= mC[15]+ d0;
2053 //fC[16]+= mC[16]+ d1;
2054 //fC[17]+= mC[17]+ d2;
2055 fC[18]+= mC[18]+ d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
2056 fC[19]+= mC[19]+ d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
2057 fC[20]+= mC[20]+ d0*mB[2][0] + d1*mB[2][1] + d2*mB[2][2];
2058
2059 d0= mB[3][0]*mV[0] + mB[3][1]*mV[1] + mB[3][2]*mV[3] - mC[21];
2060 d1= mB[3][0]*mV[1] + mB[3][1]*mV[2] + mB[3][2]*mV[4] - mC[22];
2061 d2= mB[3][0]*mV[3] + mB[3][1]*mV[4] + mB[3][2]*mV[5] - mC[23];
2062
2063 //fC[21]+= mC[21] + d0;
2064 //fC[22]+= mC[22] + d1;
2065 //fC[23]+= mC[23] + d2;
2066 fC[24]+= mC[24] + d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
2067 fC[25]+= mC[25] + d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
2068 fC[26]+= mC[26] + d0*mB[2][0] + d1*mB[2][1] + d2*mB[2][2];
2069 fC[27]+= mC[27] + d0*mB[3][0] + d1*mB[3][1] + d2*mB[3][2];
2070 }
2071
55ac3e1e 2072 SetMassConstraint(0,0);
a65041d0 2073}
2074
2075Bool_t AliKFParticleBase::InvertSym3( const Double_t A[], Double_t Ai[] )
2076{
2077 //* Invert symmetric matric stored in low-triagonal form
2078
2079 bool ret = 0;
2080 double a0 = A[0], a1 = A[1], a2 = A[2], a3 = A[3];
2081
2082 Ai[0] = a2*A[5] - A[4]*A[4];
2083 Ai[1] = a3*A[4] - a1*A[5];
2084 Ai[3] = a1*A[4] - a2*a3;
2085 Double_t det = (a0*Ai[0] + a1*Ai[1] + a3*Ai[3]);
2086 if( det>1.e-20 ) det = 1./det;
2087 else{
2088 det = 0;
2089 ret = 1;
2090 }
2091 Ai[0] *= det;
2092 Ai[1] *= det;
2093 Ai[3] *= det;
2094 Ai[2] = ( a0*A[5] - a3*a3 )*det;
2095 Ai[4] = ( a1*a3 - a0*A[4] )*det;
2096 Ai[5] = ( a0*a2 - a1*a1 )*det;
2097 return ret;
2098}
2099
f826d409 2100void AliKFParticleBase::MultQSQt( const Double_t Q[], const Double_t S[], Double_t SOut[] )
2101{
2102 //* Matrix multiplication Q*S*Q^T, Q - square matrix, S - symmetric
2103
2104 const Int_t kN= 8;
2105 Double_t mA[kN*kN];
2106
2107 for( Int_t i=0, ij=0; i<kN; i++ ){
2108 for( Int_t j=0; j<kN; j++, ++ij ){
2109 mA[ij] = 0 ;
2110 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];
2111 }
2112 }
2113
2114 for( Int_t i=0; i<kN; i++ ){
2115 for( Int_t j=0; j<=i; j++ ){
2116 Int_t ij = ( j<=i ) ? i*(i+1)/2+j :j*(j+1)/2+i;
2117 SOut[ij] = 0 ;
2118 for( Int_t k=0; k<kN; k++ ) SOut[ij] += Q[ i*kN+k ] * mA[ k*kN+j ];
2119 }
2120 }
2121}
2122
2123
2124// 72-charachters line to define the printer border
2125//3456789012345678901234567890123456789012345678901234567890123456789012
446ce366 2126