AliRoot Core  ee782a0 (ee782a0)
AliKFParticleBase.cxx
Go to the documentation of this file.
1 //---------------------------------------------------------------------------------
2 // Implementation of the AliKFParticleBase class
3 // .
4 // @author S.Gorbunov, I.Kisel, I.Kulakov, M.Zyzak
5 // @version 1.0
6 // @since 13.05.07
7 //
8 // Class to reconstruct and store the decayed particle parameters.
9 // The method is described in CBM-SOFT note 2007-003,
10 // ``Reconstruction of decayed particles based on the Kalman filter'',
11 // http://www.gsi.de/documents/DOC-2007-May-14-1.pdf
12 //
13 // This class describes general mathematics which is used by AliKFParticle class
14 //
15 // -= Copyright &copy ALICE HLT Group =-
16 //_________________________________________________________________________________
17 
18 
19 #include "AliKFParticleBase.h"
20 #include "TMath.h"
21 
22 #include <iostream>
23 ClassImp(AliKFParticleBase)
24 
25 
26 AliKFParticleBase::AliKFParticleBase() :fQ(0), fNDF(-3), fChi2(0), fSFromDecay(0), fAtProductionVertex(0), fIsLinearized(0),
27  fConstructMethod(2), SumDaughterMass(0), fMassHypo(-1)
28 {
29  //* Constructor
30 
31  Initialize();
32 }
33 
34 void AliKFParticleBase::Initialize( const Double_t Param[], const Double_t Cov[], Int_t Charge, Double_t Mass )
35 {
36  // Constructor from "cartesian" track, particle mass hypothesis should be provided
37  //
38  // Param[6] = { X, Y, Z, Px, Py, Pz } - position and momentum
39  // Cov [21] = lower-triangular part of the covariance matrix:
40  //
41  // ( 0 . . . . . )
42  // ( 1 2 . . . . )
43  // Cov. matrix = ( 3 4 5 . . . ) - numbering of covariance elements in Cov[]
44  // ( 6 7 8 9 . . )
45  // ( 10 11 12 13 14 . )
46  // ( 15 16 17 18 19 20 )
47 
48 
49  for( Int_t i=0; i<6 ; i++ ) fP[i] = Param[i];
50  for( Int_t i=0; i<21; i++ ) fC[i] = Cov[i];
51 
52  Double_t energy = TMath::Sqrt( Mass*Mass + fP[3]*fP[3] + fP[4]*fP[4] + fP[5]*fP[5]);
53  fP[6] = energy;
54  fP[7] = 0;
55  fQ = Charge;
56  fNDF = 0;
57  fChi2 = 0;
59  fIsLinearized = 0;
60  fSFromDecay = 0;
61 
62  Double_t energyInv = 1./energy;
63  Double_t
64  h0 = fP[3]*energyInv,
65  h1 = fP[4]*energyInv,
66  h2 = fP[5]*energyInv;
67 
68  fC[21] = h0*fC[ 6] + h1*fC[10] + h2*fC[15];
69  fC[22] = h0*fC[ 7] + h1*fC[11] + h2*fC[16];
70  fC[23] = h0*fC[ 8] + h1*fC[12] + h2*fC[17];
71  fC[24] = h0*fC[ 9] + h1*fC[13] + h2*fC[18];
72  fC[25] = h0*fC[13] + h1*fC[14] + h2*fC[19];
73  fC[26] = h0*fC[18] + h1*fC[19] + h2*fC[20];
74  fC[27] = ( h0*h0*fC[ 9] + h1*h1*fC[14] + h2*h2*fC[20]
75  + 2*(h0*h1*fC[13] + h0*h2*fC[18] + h1*h2*fC[19] ) );
76  for( Int_t i=28; i<36; i++ ) fC[i] = 0;
77  fC[35] = 1.;
78 
79  SumDaughterMass = Mass;
80  fMassHypo = Mass;
81 }
82 
84 {
85  //* Initialise covariance matrix and set current parameters to 0.0
86 
87  for( Int_t i=0; i<8; i++) fP[i] = 0;
88  for(Int_t i=0;i<36;++i) fC[i]=0.;
89  fC[0] = fC[2] = fC[5] = 100.;
90  fC[35] = 1.;
91  fNDF = -3;
92  fChi2 = 0.;
93  fQ = 0;
94  fSFromDecay = 0;
96  fVtxGuess[0]=fVtxGuess[1]=fVtxGuess[2]=0.;
97  fIsLinearized = 0;
98  SumDaughterMass = 0;
99  fMassHypo = -1;
100 }
101 
102 void AliKFParticleBase::SetVtxGuess( Double_t x, Double_t y, Double_t z )
103 {
104  //* Set decay vertex parameters for linearisation
105 
106  fVtxGuess[0] = x;
107  fVtxGuess[1] = y;
108  fVtxGuess[2] = z;
109  fIsLinearized = 1;
110 }
111 
112 Int_t AliKFParticleBase::GetMomentum( Double_t &p, Double_t &error ) const
113 {
114  //* Calculate particle momentum
115 
116  Double_t x = fP[3];
117  Double_t y = fP[4];
118  Double_t z = fP[5];
119 
120  Double_t x2 = x*x;
121  Double_t y2 = y*y;
122  Double_t z2 = z*z;
123  Double_t p2 = x2+y2+z2;
124  p = TMath::Sqrt(p2);
125 
126  error = (x2*fC[9]+y2*fC[14]+z2*fC[20] + 2*(x*y*fC[13]+x*z*fC[18]+y*z*fC[19]) );
127  if( error>1.e-16 && p>1.e-4 ){
128  error = TMath::Sqrt(error)/p;
129  return 0;
130  }
131  error = 1.e8;
132  return 1;
133 }
134 
135 Int_t AliKFParticleBase::GetPt( Double_t &pt, Double_t &error ) const
136 {
137  //* Calculate particle transverse momentum
138 
139  Double_t px = fP[3];
140  Double_t py = fP[4];
141  Double_t px2 = px*px;
142  Double_t py2 = py*py;
143  Double_t pt2 = px2+py2;
144  pt = TMath::Sqrt(pt2);
145  error = (px2*fC[9] + py2*fC[14] + 2*px*py*fC[13] );
146  if( error>0 && pt>1.e-4 ){
147  error = TMath::Sqrt(error)/pt;
148  return 0;
149  }
150  error = 1.e10;
151  return 1;
152 }
153 
154 Int_t AliKFParticleBase::GetEta( Double_t &eta, Double_t &error ) const
155 {
156  //* Calculate particle pseudorapidity
157 
158  Double_t px = fP[3];
159  Double_t py = fP[4];
160  Double_t pz = fP[5];
161  Double_t pt2 = px*px + py*py;
162  Double_t p2 = pt2 + pz*pz;
163  Double_t p = TMath::Sqrt(p2);
164  Double_t a = p + pz;
165  Double_t b = p - pz;
166  eta = 1.e10;
167  if( b > 1.e-8 ){
168  Double_t c = a/b;
169  if( c>1.e-8 ) eta = 0.5*TMath::Log(c);
170  }
171  Double_t h3 = -px*pz;
172  Double_t h4 = -py*pz;
173  Double_t pt4 = pt2*pt2;
174  Double_t p2pt4 = p2*pt4;
175  error = (h3*h3*fC[9] + h4*h4*fC[14] + pt4*fC[20] + 2*( h3*(h4*fC[13] + fC[18]*pt2) + pt2*h4*fC[19] ) );
176 
177  if( error>0 && p2pt4>1.e-10 ){
178  error = TMath::Sqrt(error/p2pt4);
179  return 0;
180  }
181 
182  error = 1.e10;
183  return 1;
184 }
185 
186 Int_t AliKFParticleBase::GetPhi( Double_t &phi, Double_t &error ) const
187 {
188  //* Calculate particle polar angle
189 
190  Double_t px = fP[3];
191  Double_t py = fP[4];
192  Double_t px2 = px*px;
193  Double_t py2 = py*py;
194  Double_t pt2 = px2 + py2;
195  phi = TMath::ATan2(py,px);
196  error = (py2*fC[9] + px2*fC[14] - 2*px*py*fC[13] );
197  if( error>0 && pt2>1.e-4 ){
198  error = TMath::Sqrt(error)/pt2;
199  return 0;
200  }
201  error = 1.e10;
202  return 1;
203 }
204 
205 Int_t AliKFParticleBase::GetR( Double_t &r, Double_t &error ) const
206 {
207  //* Calculate distance to the origin
208 
209  Double_t x = fP[0];
210  Double_t y = fP[1];
211  Double_t x2 = x*x;
212  Double_t y2 = y*y;
213  r = TMath::Sqrt(x2 + y2);
214  error = (x2*fC[0] + y2*fC[2] - 2*x*y*fC[1] );
215  if( error>0 && r>1.e-4 ){
216  error = TMath::Sqrt(error)/r;
217  return 0;
218  }
219  error = 1.e10;
220  return 1;
221 }
222 
223 Int_t AliKFParticleBase::GetMass( Double_t &m, Double_t &error ) const
224 {
225  //* Calculate particle mass
226 
227  // s = sigma^2 of m2/2
228 
229  Double_t s = ( fP[3]*fP[3]*fC[9] + fP[4]*fP[4]*fC[14] + fP[5]*fP[5]*fC[20]
230  + fP[6]*fP[6]*fC[27]
231  +2*( + fP[3]*fP[4]*fC[13] + fP[5]*(fP[3]*fC[18] + fP[4]*fC[19])
232  - fP[6]*( fP[3]*fC[24] + fP[4]*fC[25] + fP[5]*fC[26] ) )
233  );
234 // Double_t m2 = TMath::Abs(fP[6]*fP[6] - fP[3]*fP[3] - fP[4]*fP[4] - fP[5]*fP[5]);
235 // m = TMath::Sqrt(m2);
236 // if( m>1.e-10 ){
237 // if( s>=0 ){
238 // error = TMath::Sqrt(s)/m;
239 // return 0;
240 // }
241 // }
242 // error = 1.e20;
243 // return 1;
244  Double_t m2 = (fP[6]*fP[6] - fP[3]*fP[3] - fP[4]*fP[4] - fP[5]*fP[5]);
245 
246  if(m2<0.)
247  {
248  error = 1.e20;
249  m = -TMath::Sqrt(-m2);
250  return 1;
251  }
252 
253  m = TMath::Sqrt(m2);
254  if( m>1.e-6 ){
255  if( s >= 0 ) {
256  error = TMath::Sqrt(s)/m;
257  return 0;
258  }
259  }
260  else {
261  error = 1.e20;
262  return 0;
263  }
264  error = 1.e20;
265 
266  return 1;
267 }
268 
269 
270 Int_t AliKFParticleBase::GetDecayLength( Double_t &l, Double_t &error ) const
271 {
272  //* Calculate particle decay length [cm]
273 
274  Double_t x = fP[3];
275  Double_t y = fP[4];
276  Double_t z = fP[5];
277  Double_t t = fP[7];
278  Double_t x2 = x*x;
279  Double_t y2 = y*y;
280  Double_t z2 = z*z;
281  Double_t p2 = x2+y2+z2;
282  l = t*TMath::Sqrt(p2);
283  if( p2>1.e-4){
284  error = p2*fC[35] + t*t/p2*(x2*fC[9]+y2*fC[14]+z2*fC[20]
285  + 2*(x*y*fC[13]+x*z*fC[18]+y*z*fC[19]) )
286  + 2*t*(x*fC[31]+y*fC[32]+z*fC[33]);
287  error = TMath::Sqrt(TMath::Abs(error));
288  return 0;
289  }
290  error = 1.e20;
291  return 1;
292 }
293 
294 Int_t AliKFParticleBase::GetDecayLengthXY( Double_t &l, Double_t &error ) const
295 {
296  //* Calculate particle decay length in XY projection [cm]
297 
298  Double_t x = fP[3];
299  Double_t y = fP[4];
300  Double_t t = fP[7];
301  Double_t x2 = x*x;
302  Double_t y2 = y*y;
303  Double_t pt2 = x2+y2;
304  l = t*TMath::Sqrt(pt2);
305  if( pt2>1.e-4){
306  error = pt2*fC[35] + t*t/pt2*(x2*fC[9]+y2*fC[14] + 2*x*y*fC[13] )
307  + 2*t*(x*fC[31]+y*fC[32]);
308  error = TMath::Sqrt(TMath::Abs(error));
309  return 0;
310  }
311  error = 1.e20;
312  return 1;
313 }
314 
315 
316 Int_t AliKFParticleBase::GetLifeTime( Double_t &tauC, Double_t &error ) const
317 {
318  //* Calculate particle decay time [s]
319 
320  Double_t m, dm;
321  GetMass( m, dm );
322  Double_t cTM = (-fP[3]*fC[31] - fP[4]*fC[32] - fP[5]*fC[33] + fP[6]*fC[34]);
323  tauC = fP[7]*m;
324  error = m*m*fC[35] + 2*fP[7]*cTM + fP[7]*fP[7]*dm*dm;
325  if( error > 0 ){
326  error = TMath::Sqrt( error );
327  return 0;
328  }
329  error = 1.e20;
330  return 1;
331 }
332 
333 
335 {
336  //* Add daughter via operator+=
337 
338  AddDaughter( Daughter );
339 }
340 
341 Double_t AliKFParticleBase::GetSCorrection( const Double_t Part[], const Double_t XYZ[] )
342 {
343  //* Get big enough correction for S error to let the particle Part be fitted to XYZ point
344 
345  Double_t d[3] = { XYZ[0]-Part[0], XYZ[1]-Part[1], XYZ[2]-Part[2] };
346  Double_t p2 = Part[3]*Part[3]+Part[4]*Part[4]+Part[5]*Part[5];
347  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.;
348  return sigmaS;
349 }
350 
351 void AliKFParticleBase::GetMeasurement( const Double_t XYZ[], Double_t m[], Double_t V[] ) const
352 {
353  //* Get additional covariances V used during measurement
354 
355  Double_t b[3];
356  GetFieldValue( XYZ, b );
357  const Double_t kCLight = 0.000299792458;
358  b[0]*=kCLight; b[1]*=kCLight; b[2]*=kCLight;
359 
360  Transport( GetDStoPoint(XYZ), m, V );
361 
362  Double_t sigmaS = GetSCorrection( m, XYZ );
363 
364  Double_t h[6];
365 
366  h[0] = m[3]*sigmaS;
367  h[1] = m[4]*sigmaS;
368  h[2] = m[5]*sigmaS;
369  h[3] = ( h[1]*b[2]-h[2]*b[1] )*GetQ();
370  h[4] = ( h[2]*b[0]-h[0]*b[2] )*GetQ();
371  h[5] = ( h[0]*b[1]-h[1]*b[0] )*GetQ();
372 
373  V[ 0]+= h[0]*h[0];
374  V[ 1]+= h[1]*h[0];
375  V[ 2]+= h[1]*h[1];
376  V[ 3]+= h[2]*h[0];
377  V[ 4]+= h[2]*h[1];
378  V[ 5]+= h[2]*h[2];
379 
380  V[ 6]+= h[3]*h[0];
381  V[ 7]+= h[3]*h[1];
382  V[ 8]+= h[3]*h[2];
383  V[ 9]+= h[3]*h[3];
384 
385  V[10]+= h[4]*h[0];
386  V[11]+= h[4]*h[1];
387  V[12]+= h[4]*h[2];
388  V[13]+= h[4]*h[3];
389  V[14]+= h[4]*h[4];
390 
391  V[15]+= h[5]*h[0];
392  V[16]+= h[5]*h[1];
393  V[17]+= h[5]*h[2];
394  V[18]+= h[5]*h[3];
395  V[19]+= h[5]*h[4];
396  V[20]+= h[5]*h[5];
397 }
398 
400 {
401  if( fNDF<-1 ){ // first daughter -> just copy
402  fNDF = -1;
403  fQ = Daughter.GetQ();
404  for( Int_t i=0; i<7; i++) fP[i] = Daughter.fP[i];
405  for( Int_t i=0; i<28; i++) fC[i] = Daughter.fC[i];
406  fSFromDecay = 0;
407  fMassHypo = Daughter.fMassHypo;
408  SumDaughterMass = Daughter.SumDaughterMass;
409  return;
410  }
411 
412  if(fConstructMethod == 0)
413  AddDaughterWithEnergyFit(Daughter);
414  else if(fConstructMethod == 1)
415  AddDaughterWithEnergyCalc(Daughter);
416  else if(fConstructMethod == 2)
417  AddDaughterWithEnergyFitMC(Daughter);
418 
419  SumDaughterMass += Daughter.SumDaughterMass;
420  fMassHypo = -1;
421 }
422 
424 {
425  //* Energy considered as an independent veriable, fitted independently from momentum, without any constraints on mass
426 
427  //* Add daughter
428 
430 
431  Double_t b[3];
432  Int_t maxIter = 1;
433 
434  if( !fIsLinearized ){
435  if( fNDF==-1 ){
436  Double_t ds, ds1;
437  GetDStoParticle(Daughter, ds, ds1);
438  TransportToDS( ds );
439  Double_t m[8];
440  Double_t mCd[36];
441  Daughter.Transport( ds1, m, mCd );
442  fVtxGuess[0] = .5*( fP[0] + m[0] );
443  fVtxGuess[1] = .5*( fP[1] + m[1] );
444  fVtxGuess[2] = .5*( fP[2] + m[2] );
445  } else {
446  fVtxGuess[0] = fP[0];
447  fVtxGuess[1] = fP[1];
448  fVtxGuess[2] = fP[2];
449  }
450  maxIter = 3;
451  }
452 
453  for( Int_t iter=0; iter<maxIter; iter++ ){
454 
455  {
456  GetFieldValue( fVtxGuess, b );
457  const Double_t kCLight = 0.000299792458;
458  b[0]*=kCLight; b[1]*=kCLight; b[2]*=kCLight;
459  }
460 
461  Double_t *ffP = fP, *ffC = fC, tmpP[8], tmpC[36];
462  if( fNDF==-1 ){
463  GetMeasurement( fVtxGuess, tmpP, tmpC );
464  ffP = tmpP;
465  ffC = tmpC;
466  }
467 
468  Double_t m[8], mV[36];
469 
470  if( Daughter.fC[35]>0 ){
471  Daughter.GetMeasurement( fVtxGuess, m, mV );
472  } else {
473  for( Int_t i=0; i<8; i++ ) m[i] = Daughter.fP[i];
474  for( Int_t i=0; i<36; i++ ) mV[i] = Daughter.fC[i];
475  }
476  //*
477 
478  Double_t mS[6];
479  {
480  Double_t mSi[6] = { ffC[0]+mV[0],
481  ffC[1]+mV[1], ffC[2]+mV[2],
482  ffC[3]+mV[3], ffC[4]+mV[4], ffC[5]+mV[5] };
483 
484  mS[0] = mSi[2]*mSi[5] - mSi[4]*mSi[4];
485  mS[1] = mSi[3]*mSi[4] - mSi[1]*mSi[5];
486  mS[2] = mSi[0]*mSi[5] - mSi[3]*mSi[3];
487  mS[3] = mSi[1]*mSi[4] - mSi[2]*mSi[3];
488  mS[4] = mSi[1]*mSi[3] - mSi[0]*mSi[4];
489  mS[5] = mSi[0]*mSi[2] - mSi[1]*mSi[1];
490 
491  Double_t s = ( mSi[0]*mS[0] + mSi[1]*mS[1] + mSi[3]*mS[3] );
492  s = ( TMath::Abs(s) > 1.E-20 ) ?1./s :0;
493  mS[0]*=s;
494  mS[1]*=s;
495  mS[2]*=s;
496  mS[3]*=s;
497  mS[4]*=s;
498  mS[5]*=s;
499  }
500  //* Residual (measured - estimated)
501 
502  Double_t zeta[3] = { m[0]-ffP[0], m[1]-ffP[1], m[2]-ffP[2] };
503 
504  //* CHt = CH' - D'
505 
506  Double_t mCHt0[7], mCHt1[7], mCHt2[7];
507 
508  mCHt0[0]=ffC[ 0] ; mCHt1[0]=ffC[ 1] ; mCHt2[0]=ffC[ 3] ;
509  mCHt0[1]=ffC[ 1] ; mCHt1[1]=ffC[ 2] ; mCHt2[1]=ffC[ 4] ;
510  mCHt0[2]=ffC[ 3] ; mCHt1[2]=ffC[ 4] ; mCHt2[2]=ffC[ 5] ;
511  mCHt0[3]=ffC[ 6]-mV[ 6]; mCHt1[3]=ffC[ 7]-mV[ 7]; mCHt2[3]=ffC[ 8]-mV[ 8];
512  mCHt0[4]=ffC[10]-mV[10]; mCHt1[4]=ffC[11]-mV[11]; mCHt2[4]=ffC[12]-mV[12];
513  mCHt0[5]=ffC[15]-mV[15]; mCHt1[5]=ffC[16]-mV[16]; mCHt2[5]=ffC[17]-mV[17];
514  mCHt0[6]=ffC[21]-mV[21]; mCHt1[6]=ffC[22]-mV[22]; mCHt2[6]=ffC[23]-mV[23];
515 
516  //* Kalman gain K = mCH'*S
517 
518  Double_t k0[7], k1[7], k2[7];
519 
520  for(Int_t i=0;i<7;++i){
521  k0[i] = mCHt0[i]*mS[0] + mCHt1[i]*mS[1] + mCHt2[i]*mS[3];
522  k1[i] = mCHt0[i]*mS[1] + mCHt1[i]*mS[2] + mCHt2[i]*mS[4];
523  k2[i] = mCHt0[i]*mS[3] + mCHt1[i]*mS[4] + mCHt2[i]*mS[5];
524  }
525 
526  //* New estimation of the vertex position
527 
528  if( iter<maxIter-1 ){
529  for(Int_t i=0; i<3; ++i)
530  fVtxGuess[i]= ffP[i] + k0[i]*zeta[0]+k1[i]*zeta[1]+k2[i]*zeta[2];
531  continue;
532  }
533 
534  // last itearation -> update the particle
535 
536  //* Add the daughter momentum to the particle momentum
537 
538  ffP[ 3] += m[ 3];
539  ffP[ 4] += m[ 4];
540  ffP[ 5] += m[ 5];
541  ffP[ 6] += m[ 6];
542 
543  ffC[ 9] += mV[ 9];
544  ffC[13] += mV[13];
545  ffC[14] += mV[14];
546  ffC[18] += mV[18];
547  ffC[19] += mV[19];
548  ffC[20] += mV[20];
549  ffC[24] += mV[24];
550  ffC[25] += mV[25];
551  ffC[26] += mV[26];
552  ffC[27] += mV[27];
553 
554 
555  //* New estimation of the vertex position r += K*zeta
556 
557  for(Int_t i=0;i<7;++i)
558  fP[i] = ffP[i] + k0[i]*zeta[0] + k1[i]*zeta[1] + k2[i]*zeta[2];
559 
560  //* New covariance matrix C -= K*(mCH')'
561 
562  for(Int_t i=0, k=0;i<7;++i){
563  for(Int_t j=0;j<=i;++j,++k){
564  fC[k] = ffC[k] - (k0[i]*mCHt0[j] + k1[i]*mCHt1[j] + k2[i]*mCHt2[j] );
565  }
566  }
567 
568  //* Calculate Chi^2
569 
570  fNDF += 2;
571  fQ += Daughter.GetQ();
572  fSFromDecay = 0;
573  fChi2 += (mS[0]*zeta[0] + mS[1]*zeta[1] + mS[3]*zeta[2])*zeta[0]
574  + (mS[1]*zeta[0] + mS[2]*zeta[1] + mS[4]*zeta[2])*zeta[1]
575  + (mS[3]*zeta[0] + mS[4]*zeta[1] + mS[5]*zeta[2])*zeta[2];
576 
577  }
578 }
579 
581 {
582  //* Energy considered as a dependent variable, calculated from the momentum and mass hypothesis
583 
584  //* Add daughter
585 
587 
588  Double_t b[3];
589  Int_t maxIter = 1;
590 
591  if( !fIsLinearized ){
592  if( fNDF==-1 ){
593  Double_t ds, ds1;
594  GetDStoParticle(Daughter, ds, ds1);
595  TransportToDS( ds );
596  Double_t m[8];
597  Double_t mCd[36];
598  Daughter.Transport( ds1, m, mCd );
599  fVtxGuess[0] = .5*( fP[0] + m[0] );
600  fVtxGuess[1] = .5*( fP[1] + m[1] );
601  fVtxGuess[2] = .5*( fP[2] + m[2] );
602  } else {
603  fVtxGuess[0] = fP[0];
604  fVtxGuess[1] = fP[1];
605  fVtxGuess[2] = fP[2];
606  }
607  maxIter = 3;
608  }
609 
610  for( Int_t iter=0; iter<maxIter; iter++ ){
611 
612  {
613  GetFieldValue( fVtxGuess, b );
614  const Double_t kCLight = 0.000299792458;
615  b[0]*=kCLight; b[1]*=kCLight; b[2]*=kCLight;
616  }
617 
618  Double_t *ffP = fP, *ffC = fC, tmpP[8], tmpC[36];
619  if( fNDF==-1 ){
620  GetMeasurement( fVtxGuess, tmpP, tmpC );
621  ffP = tmpP;
622  ffC = tmpC;
623  }
624 
625  Double_t m[8], mV[36];
626 
627  if( Daughter.fC[35]>0 ){
628  Daughter.GetMeasurement( fVtxGuess, m, mV );
629  } else {
630  for( Int_t i=0; i<8; i++ ) m[i] = Daughter.fP[i];
631  for( Int_t i=0; i<36; i++ ) mV[i] = Daughter.fC[i];
632  }
633 
634  double massMf2 = m[6]*m[6] - (m[3]*m[3] + m[4]*m[4] + m[5]*m[5]);
635  double massRf2 = fP[6]*fP[6] - (fP[3]*fP[3] + fP[4]*fP[4] + fP[5]*fP[5]);
636 
637  //*
638 
639  Double_t mS[6];
640  {
641  Double_t mSi[6] = { ffC[0]+mV[0],
642  ffC[1]+mV[1], ffC[2]+mV[2],
643  ffC[3]+mV[3], ffC[4]+mV[4], ffC[5]+mV[5] };
644 
645  mS[0] = mSi[2]*mSi[5] - mSi[4]*mSi[4];
646  mS[1] = mSi[3]*mSi[4] - mSi[1]*mSi[5];
647  mS[2] = mSi[0]*mSi[5] - mSi[3]*mSi[3];
648  mS[3] = mSi[1]*mSi[4] - mSi[2]*mSi[3];
649  mS[4] = mSi[1]*mSi[3] - mSi[0]*mSi[4];
650  mS[5] = mSi[0]*mSi[2] - mSi[1]*mSi[1];
651 
652  Double_t s = ( mSi[0]*mS[0] + mSi[1]*mS[1] + mSi[3]*mS[3] );
653 
654  s = ( s > 1.E-20 ) ?1./s :0;
655  mS[0]*=s;
656  mS[1]*=s;
657  mS[2]*=s;
658  mS[3]*=s;
659  mS[4]*=s;
660  mS[5]*=s;
661  }
662 
663  //* Residual (measured - estimated)
664 
665  Double_t zeta[3] = { m[0]-ffP[0], m[1]-ffP[1], m[2]-ffP[2] };
666 
667  //* CHt = CH' - D'
668 
669  Double_t mCHt0[6], mCHt1[6], mCHt2[6];
670 
671  mCHt0[0]=ffC[ 0] ; mCHt1[0]=ffC[ 1] ; mCHt2[0]=ffC[ 3] ;
672  mCHt0[1]=ffC[ 1] ; mCHt1[1]=ffC[ 2] ; mCHt2[1]=ffC[ 4] ;
673  mCHt0[2]=ffC[ 3] ; mCHt1[2]=ffC[ 4] ; mCHt2[2]=ffC[ 5] ;
674  mCHt0[3]=ffC[ 6]-mV[ 6]; mCHt1[3]=ffC[ 7]-mV[ 7]; mCHt2[3]=ffC[ 8]-mV[ 8];
675  mCHt0[4]=ffC[10]-mV[10]; mCHt1[4]=ffC[11]-mV[11]; mCHt2[4]=ffC[12]-mV[12];
676  mCHt0[5]=ffC[15]-mV[15]; mCHt1[5]=ffC[16]-mV[16]; mCHt2[5]=ffC[17]-mV[17];
677 
678  //* Kalman gain K = mCH'*S
679 
680  Double_t k0[6], k1[6], k2[6];
681 
682  for(Int_t i=0;i<6;++i){
683  k0[i] = mCHt0[i]*mS[0] + mCHt1[i]*mS[1] + mCHt2[i]*mS[3];
684  k1[i] = mCHt0[i]*mS[1] + mCHt1[i]*mS[2] + mCHt2[i]*mS[4];
685  k2[i] = mCHt0[i]*mS[3] + mCHt1[i]*mS[4] + mCHt2[i]*mS[5];
686  }
687 
688  //* New estimation of the vertex position
689 
690  if( iter<maxIter-1 ){
691  for(Int_t i=0; i<3; ++i)
692  fVtxGuess[i]= ffP[i] + k0[i]*zeta[0]+k1[i]*zeta[1]+k2[i]*zeta[2];
693  continue;
694  }
695 
696  //* find mf and mVf - optimum value of the measurement and its covariance matrix
697  //* mVHt = V*H'
698  Double_t mVHt0[6], mVHt1[6], mVHt2[6];
699 
700  mVHt0[0]= mV[ 0] ; mVHt1[0]= mV[ 1] ; mVHt2[0]= mV[ 3] ;
701  mVHt0[1]= mV[ 1] ; mVHt1[1]= mV[ 2] ; mVHt2[1]= mV[ 4] ;
702  mVHt0[2]= mV[ 3] ; mVHt1[2]= mV[ 4] ; mVHt2[2]= mV[ 5] ;
703  mVHt0[3]= mV[ 6] ; mVHt1[3]= mV[ 7] ; mVHt2[3]= mV[ 8] ;
704  mVHt0[4]= mV[10] ; mVHt1[4]= mV[11] ; mVHt2[4]= mV[12] ;
705  mVHt0[5]= mV[15] ; mVHt1[5]= mV[16] ; mVHt2[5]= mV[17] ;
706 
707  //* Kalman gain Km = mCH'*S
708 
709  Double_t km0[6], km1[6], km2[6];
710 
711  for(Int_t i=0;i<6;++i){
712  km0[i] = mVHt0[i]*mS[0] + mVHt1[i]*mS[1] + mVHt2[i]*mS[3];
713  km1[i] = mVHt0[i]*mS[1] + mVHt1[i]*mS[2] + mVHt2[i]*mS[4];
714  km2[i] = mVHt0[i]*mS[3] + mVHt1[i]*mS[4] + mVHt2[i]*mS[5];
715  }
716 
717  Double_t mf[7] = { m[0], m[1], m[2], m[3], m[4], m[5], m[6] };
718 
719  for(Int_t i=0;i<6;++i)
720  mf[i] = mf[i] - km0[i]*zeta[0] - km1[i]*zeta[1] - km2[i]*zeta[2];
721 
722  Double_t energyMf = TMath::Sqrt( massMf2 + (mf[3]*mf[3] + mf[4]*mf[4] + mf[5]*mf[5]) );
723 
724  Double_t mVf[28];
725  for(Int_t iC=0; iC<28; iC++)
726  mVf[iC] = mV[iC];
727 
728  //* hmf = d(energyMf)/d(mf)
729  Double_t hmf[7];
730  if( TMath::Abs(energyMf) < 1.e-10) hmf[3] = 0; else hmf[3] = mf[3]/energyMf;
731  if( TMath::Abs(energyMf) < 1.e-10) hmf[4] = 0; else hmf[4] = mf[4]/energyMf;
732  if( TMath::Abs(energyMf) < 1.e-10) hmf[5] = 0; else hmf[5] = mf[5]/energyMf;
733 // if( TMath::Abs(energyMf) < 1.e-10) hmf[6] = 0; else hmf[6] = mf[6]/energyMf;
734  hmf[6] = 0;
735 
736  for(Int_t i=0, k=0;i<6;++i){
737  for(Int_t j=0;j<=i;++j,++k){
738  mVf[k] = mVf[k] - (km0[i]*mVHt0[j] + km1[i]*mVHt1[j] + km2[i]*mVHt2[j] );
739  }
740  }
741  Double_t mVf24 = mVf[24], mVf25 = mVf[25], mVf26 = mVf[26];
742  mVf[21] = mVf[6 ]*hmf[3] + mVf[10]*hmf[4] + mVf[15]*hmf[5] + mVf[21]*hmf[6];
743  mVf[22] = mVf[7 ]*hmf[3] + mVf[11]*hmf[4] + mVf[16]*hmf[5] + mVf[22]*hmf[6];
744  mVf[23] = mVf[8 ]*hmf[3] + mVf[12]*hmf[4] + mVf[17]*hmf[5] + mVf[23]*hmf[6];
745  mVf[24] = mVf[9 ]*hmf[3] + mVf[13]*hmf[4] + mVf[18]*hmf[5] + mVf[24]*hmf[6];
746  mVf[25] = mVf[13]*hmf[3] + mVf[14]*hmf[4] + mVf[19]*hmf[5] + mVf[25]*hmf[6];
747  mVf[26] = mVf[18]*hmf[3] + mVf[19]*hmf[4] + mVf[20]*hmf[5] + mVf[26]*hmf[6];
748  mVf[27] = mVf[24]*hmf[3] + mVf[25]*hmf[4] + mVf[26]*hmf[5] + (mVf24*hmf[3] + mVf25*hmf[4] + mVf26*hmf[5] + mVf[27]*hmf[6])*hmf[6]; //here mVf[] are already modified
749 
750  mf[6] = energyMf;
751 
752  //* find rf and mCf - optimum value of the measurement and its covariance matrix
753 
754  //* mCCHt = C*H'
755  Double_t mCCHt0[6], mCCHt1[6], mCCHt2[6];
756 
757  mCCHt0[0]=ffC[ 0]; mCCHt1[0]=ffC[ 1]; mCCHt2[0]=ffC[ 3];
758  mCCHt0[1]=ffC[ 1]; mCCHt1[1]=ffC[ 2]; mCCHt2[1]=ffC[ 4];
759  mCCHt0[2]=ffC[ 3]; mCCHt1[2]=ffC[ 4]; mCCHt2[2]=ffC[ 5];
760  mCCHt0[3]=ffC[ 6]; mCCHt1[3]=ffC[ 7]; mCCHt2[3]=ffC[ 8];
761  mCCHt0[4]=ffC[10]; mCCHt1[4]=ffC[11]; mCCHt2[4]=ffC[12];
762  mCCHt0[5]=ffC[15]; mCCHt1[5]=ffC[16]; mCCHt2[5]=ffC[17];
763 
764  //* Kalman gain Krf = mCH'*S
765 
766  Double_t krf0[6], krf1[6], krf2[6];
767 
768  for(Int_t i=0;i<6;++i){
769  krf0[i] = mCCHt0[i]*mS[0] + mCCHt1[i]*mS[1] + mCCHt2[i]*mS[3];
770  krf1[i] = mCCHt0[i]*mS[1] + mCCHt1[i]*mS[2] + mCCHt2[i]*mS[4];
771  krf2[i] = mCCHt0[i]*mS[3] + mCCHt1[i]*mS[4] + mCCHt2[i]*mS[5];
772  }
773  Double_t rf[7] = { ffP[0], ffP[1], ffP[2], ffP[3], ffP[4], ffP[5], ffP[6] };
774 
775  for(Int_t i=0;i<6;++i)
776  rf[i] = rf[i] + krf0[i]*zeta[0] + krf1[i]*zeta[1] + krf2[i]*zeta[2];
777 
778  Double_t energyRf = TMath::Sqrt( massRf2 + (rf[3]*rf[3] + rf[4]*rf[4] + rf[5]*rf[5]) );
779 
780  Double_t mCf[28];
781  for(Int_t iC=0; iC<28; iC++)
782  mCf[iC] = ffC[iC];
783  //* hrf = d(Erf)/d(rf)
784  Double_t hrf[7];
785  if( TMath::Abs(energyRf) < 1.e-10) hrf[3] = 0; else hrf[3] = rf[3]/energyRf;
786  if( TMath::Abs(energyRf) < 1.e-10) hrf[4] = 0; else hrf[4] = rf[4]/energyRf;
787  if( TMath::Abs(energyRf) < 1.e-10) hrf[5] = 0; else hrf[5] = rf[5]/energyRf;
788 // if( TMath::Abs(energyRf) < 1.e-10) hrf[6] = 0; else hrf[6] = rf[6]/energyRf;
789  hrf[6] = 0;
790 
791  for(Int_t i=0, k=0;i<6;++i){
792  for(Int_t j=0;j<=i;++j,++k){
793  mCf[k] = mCf[k] - (krf0[i]*mCCHt0[j] + krf1[i]*mCCHt1[j] + krf2[i]*mCCHt2[j] );
794  }
795  }
796  Double_t mCf24 = mCf[24], mCf25 = mCf[25], mCf26 = mCf[26];
797  mCf[21] = mCf[6 ]*hrf[3] + mCf[10]*hrf[4] + mCf[15]*hrf[5] + mCf[21]*hrf[6];
798  mCf[22] = mCf[7 ]*hrf[3] + mCf[11]*hrf[4] + mCf[16]*hrf[5] + mCf[22]*hrf[6];
799  mCf[23] = mCf[8 ]*hrf[3] + mCf[12]*hrf[4] + mCf[17]*hrf[5] + mCf[23]*hrf[6];
800  mCf[24] = mCf[9 ]*hrf[3] + mCf[13]*hrf[4] + mCf[18]*hrf[5] + mCf[24]*hrf[6];
801  mCf[25] = mCf[13]*hrf[3] + mCf[14]*hrf[4] + mCf[19]*hrf[5] + mCf[25]*hrf[6];
802  mCf[26] = mCf[18]*hrf[3] + mCf[19]*hrf[4] + mCf[20]*hrf[5] + mCf[26]*hrf[6];
803  mCf[27] = mCf[24]*hrf[3] + mCf[25]*hrf[4] + mCf[26]*hrf[5] + (mCf24*hrf[3] + mCf25*hrf[4] + mCf26*hrf[5] + mCf[27]*hrf[6])*hrf[6]; //here mCf[] are already modified
804 
805  for(Int_t iC=21; iC<28; iC++)
806  {
807  ffC[iC] = mCf[iC];
808  mV[iC] = mVf[iC];
809  }
810 
811  fP[6] = energyRf + energyMf;
812  rf[6] = energyRf;
813 
814  //Double_t Dvv[3][3]; do not need this
815  Double_t mDvp[3][3];
816  // Double_t mDpv[3][3];
817  Double_t mDpp[3][3];
818  Double_t mDe[7];
819 
820  for(int i=0; i<3; i++)
821  {
822  for(int j=0; j<3; j++)
823  {
824  mDvp[i][j] = km0[i+3]*mCCHt0[j] + km1[i+3]*mCCHt1[j] + km2[i+3]*mCCHt2[j];
825  // mDpv[i][j] = km0[i]*mCCHt0[j+3] + km1[i]*mCCHt1[j+3] + km2[i]*mCCHt2[j+3];
826  mDpp[i][j] = km0[i+3]*mCCHt0[j+3] + km1[i+3]*mCCHt1[j+3] + km2[i+3]*mCCHt2[j+3];
827  }
828  }
829 
830  mDe[0] = hmf[3]*mDvp[0][0] + hmf[4]*mDvp[1][0] + hmf[5]*mDvp[2][0];
831  mDe[1] = hmf[3]*mDvp[0][1] + hmf[4]*mDvp[1][1] + hmf[5]*mDvp[2][1];
832  mDe[2] = hmf[3]*mDvp[0][2] + hmf[4]*mDvp[1][2] + hmf[5]*mDvp[2][2];
833  mDe[3] = hmf[3]*mDpp[0][0] + hmf[4]*mDpp[1][0] + hmf[5]*mDpp[2][0];
834  mDe[4] = hmf[3]*mDpp[0][1] + hmf[4]*mDpp[1][1] + hmf[5]*mDpp[2][1];
835  mDe[5] = hmf[3]*mDpp[0][2] + hmf[4]*mDpp[1][2] + hmf[5]*mDpp[2][2];
836  mDe[6] = 2*(mDe[3]*hrf[3] + mDe[4]*hrf[4] + mDe[5]*hrf[5]);
837 
838  // last itearation -> update the particle
839 
840  //* Add the daughter momentum to the particle momentum
841 
842  ffP[ 3] += m[ 3];
843  ffP[ 4] += m[ 4];
844  ffP[ 5] += m[ 5];
845 
846  ffC[ 9] += mV[ 9];
847  ffC[13] += mV[13];
848  ffC[14] += mV[14];
849  ffC[18] += mV[18];
850  ffC[19] += mV[19];
851  ffC[20] += mV[20];
852  ffC[24] += mV[24];
853  ffC[25] += mV[25];
854  ffC[26] += mV[26];
855  ffC[27] += mV[27];
856 
857  ffC[21] += mDe[0];
858  ffC[22] += mDe[1];
859  ffC[23] += mDe[2];
860  ffC[24] += mDe[3];
861  ffC[25] += mDe[4];
862  ffC[26] += mDe[5];
863  ffC[27] += mDe[6];
864 
865  //* New estimation of the vertex position r += K*zeta
866 
867  for(Int_t i=0;i<6;++i)
868  fP[i] = ffP[i] + k0[i]*zeta[0] + k1[i]*zeta[1] + k2[i]*zeta[2];
869 
870  //* New covariance matrix C -= K*(mCH')'
871 
872  for(Int_t i=0, k=0;i<6;++i){
873  for(Int_t j=0;j<=i;++j,++k){
874  fC[k] = ffC[k] - (k0[i]*mCHt0[j] + k1[i]*mCHt1[j] + k2[i]*mCHt2[j] );
875  }
876  }
877 
878  for(int i=21; i<28; i++) fC[i] = ffC[i];
879 
880  //* Calculate Chi^2
881 
882  fNDF += 2;
883  fQ += Daughter.GetQ();
884  fSFromDecay = 0;
885  fChi2 += (mS[0]*zeta[0] + mS[1]*zeta[1] + mS[3]*zeta[2])*zeta[0]
886  + (mS[1]*zeta[0] + mS[2]*zeta[1] + mS[4]*zeta[2])*zeta[1]
887  + (mS[3]*zeta[0] + mS[4]*zeta[1] + mS[5]*zeta[2])*zeta[2];
888  }
889 }
890 
892 {
893  //* Energy considered as an independent variable, fitted independently from momentum, without any constraints on mass
894 
895  //* Add daughter
896 
898 
899  Double_t b[3];
900  Int_t maxIter = 1;
901 
902  if( !fIsLinearized ){
903  if( fNDF==-1 ){
904  Double_t ds, ds1;
905  GetDStoParticle(Daughter, ds, ds1);
906  TransportToDS( ds );
907  Double_t m[8];
908  Double_t mCd[36];
909  Daughter.Transport( ds1, m, mCd );
910  fVtxGuess[0] = .5*( fP[0] + m[0] );
911  fVtxGuess[1] = .5*( fP[1] + m[1] );
912  fVtxGuess[2] = .5*( fP[2] + m[2] );
913  } else {
914  fVtxGuess[0] = fP[0];
915  fVtxGuess[1] = fP[1];
916  fVtxGuess[2] = fP[2];
917  }
918  maxIter = 3;
919  }
920 
921  for( Int_t iter=0; iter<maxIter; iter++ ){
922 
923  {
924  GetFieldValue( fVtxGuess, b );
925  const Double_t kCLight = 0.000299792458;
926  b[0]*=kCLight; b[1]*=kCLight; b[2]*=kCLight;
927  }
928 
929  Double_t *ffP = fP, *ffC = fC, tmpP[8], tmpC[36];
930  if( fNDF==-1 ){
931  GetMeasurement( fVtxGuess, tmpP, tmpC );
932  ffP = tmpP;
933  ffC = tmpC;
934  }
935  Double_t m[8], mV[36];
936 
937  if( Daughter.fC[35]>0 ){
938  Daughter.GetMeasurement( fVtxGuess, m, mV );
939  } else {
940  for( Int_t i=0; i<8; i++ ) m[i] = Daughter.fP[i];
941  for( Int_t i=0; i<36; i++ ) mV[i] = Daughter.fC[i];
942  }
943  //*
944 
945  Double_t mS[6];
946  {
947  Double_t mSi[6] = { ffC[0]+mV[0],
948  ffC[1]+mV[1], ffC[2]+mV[2],
949  ffC[3]+mV[3], ffC[4]+mV[4], ffC[5]+mV[5] };
950 
951  mS[0] = mSi[2]*mSi[5] - mSi[4]*mSi[4];
952  mS[1] = mSi[3]*mSi[4] - mSi[1]*mSi[5];
953  mS[2] = mSi[0]*mSi[5] - mSi[3]*mSi[3];
954  mS[3] = mSi[1]*mSi[4] - mSi[2]*mSi[3];
955  mS[4] = mSi[1]*mSi[3] - mSi[0]*mSi[4];
956  mS[5] = mSi[0]*mSi[2] - mSi[1]*mSi[1];
957 
958  Double_t s = ( mSi[0]*mS[0] + mSi[1]*mS[1] + mSi[3]*mS[3] );
959 
960  s = ( s > 1.E-20 ) ?1./s :0;
961  mS[0]*=s;
962  mS[1]*=s;
963  mS[2]*=s;
964  mS[3]*=s;
965  mS[4]*=s;
966  mS[5]*=s;
967  }
968  //* Residual (measured - estimated)
969 
970  Double_t zeta[3] = { m[0]-ffP[0], m[1]-ffP[1], m[2]-ffP[2] };
971 
972 
973  //* CHt = CH'
974 
975  Double_t mCHt0[7], mCHt1[7], mCHt2[7];
976 
977  mCHt0[0]=ffC[ 0] ; mCHt1[0]=ffC[ 1] ; mCHt2[0]=ffC[ 3] ;
978  mCHt0[1]=ffC[ 1] ; mCHt1[1]=ffC[ 2] ; mCHt2[1]=ffC[ 4] ;
979  mCHt0[2]=ffC[ 3] ; mCHt1[2]=ffC[ 4] ; mCHt2[2]=ffC[ 5] ;
980  mCHt0[3]=ffC[ 6] ; mCHt1[3]=ffC[ 7] ; mCHt2[3]=ffC[ 8] ;
981  mCHt0[4]=ffC[10] ; mCHt1[4]=ffC[11] ; mCHt2[4]=ffC[12] ;
982  mCHt0[5]=ffC[15] ; mCHt1[5]=ffC[16] ; mCHt2[5]=ffC[17] ;
983  mCHt0[6]=ffC[21] ; mCHt1[6]=ffC[22] ; mCHt2[6]=ffC[23] ;
984 
985  //* Kalman gain K = mCH'*S
986 
987  Double_t k0[7], k1[7], k2[7];
988 
989  for(Int_t i=0;i<7;++i){
990  k0[i] = mCHt0[i]*mS[0] + mCHt1[i]*mS[1] + mCHt2[i]*mS[3];
991  k1[i] = mCHt0[i]*mS[1] + mCHt1[i]*mS[2] + mCHt2[i]*mS[4];
992  k2[i] = mCHt0[i]*mS[3] + mCHt1[i]*mS[4] + mCHt2[i]*mS[5];
993  }
994 
995  //* New estimation of the vertex position
996 
997  if( iter<maxIter-1 ){
998  for(Int_t i=0; i<3; ++i)
999  fVtxGuess[i]= ffP[i] + k0[i]*zeta[0]+k1[i]*zeta[1]+k2[i]*zeta[2];
1000  continue;
1001  }
1002 
1003  // last itearation -> update the particle
1004 
1005  //* VHt = VH'
1006 
1007  Double_t mVHt0[7], mVHt1[7], mVHt2[7];
1008 
1009  mVHt0[0]=mV[ 0] ; mVHt1[0]=mV[ 1] ; mVHt2[0]=mV[ 3] ;
1010  mVHt0[1]=mV[ 1] ; mVHt1[1]=mV[ 2] ; mVHt2[1]=mV[ 4] ;
1011  mVHt0[2]=mV[ 3] ; mVHt1[2]=mV[ 4] ; mVHt2[2]=mV[ 5] ;
1012  mVHt0[3]=mV[ 6] ; mVHt1[3]=mV[ 7] ; mVHt2[3]=mV[ 8] ;
1013  mVHt0[4]=mV[10] ; mVHt1[4]=mV[11] ; mVHt2[4]=mV[12] ;
1014  mVHt0[5]=mV[15] ; mVHt1[5]=mV[16] ; mVHt2[5]=mV[17] ;
1015  mVHt0[6]=mV[21] ; mVHt1[6]=mV[22] ; mVHt2[6]=mV[23] ;
1016 
1017  //* Kalman gain Km = mCH'*S
1018 
1019  Double_t km0[7], km1[7], km2[7];
1020 
1021  for(Int_t i=0;i<7;++i){
1022  km0[i] = mVHt0[i]*mS[0] + mVHt1[i]*mS[1] + mVHt2[i]*mS[3];
1023  km1[i] = mVHt0[i]*mS[1] + mVHt1[i]*mS[2] + mVHt2[i]*mS[4];
1024  km2[i] = mVHt0[i]*mS[3] + mVHt1[i]*mS[4] + mVHt2[i]*mS[5];
1025  }
1026 
1027  for(Int_t i=0;i<7;++i)
1028  ffP[i] = ffP[i] + k0[i]*zeta[0] + k1[i]*zeta[1] + k2[i]*zeta[2];
1029 
1030  for(Int_t i=0;i<7;++i)
1031  m[i] = m[i] - km0[i]*zeta[0] - km1[i]*zeta[1] - km2[i]*zeta[2];
1032 
1033  for(Int_t i=0, k=0;i<7;++i){
1034  for(Int_t j=0;j<=i;++j,++k){
1035  ffC[k] = ffC[k] - (k0[i]*mCHt0[j] + k1[i]*mCHt1[j] + k2[i]*mCHt2[j] );
1036  }
1037  }
1038 
1039  for(Int_t i=0, k=0;i<7;++i){
1040  for(Int_t j=0;j<=i;++j,++k){
1041  mV[k] = mV[k] - (km0[i]*mVHt0[j] + km1[i]*mVHt1[j] + km2[i]*mVHt2[j] );
1042  }
1043  }
1044 
1045  Double_t mDf[7][7];
1046 
1047  for(Int_t i=0;i<7;++i){
1048  for(Int_t j=0;j<7;++j){
1049  mDf[i][j] = (km0[i]*mCHt0[j] + km1[i]*mCHt1[j] + km2[i]*mCHt2[j] );
1050  }
1051  }
1052 
1053  Double_t mJ1[7][7], mJ2[7][7];
1054  for(Int_t iPar1=0; iPar1<7; iPar1++)
1055  {
1056  for(Int_t iPar2=0; iPar2<7; iPar2++)
1057  {
1058  mJ1[iPar1][iPar2] = 0;
1059  mJ2[iPar1][iPar2] = 0;
1060  }
1061  }
1062 
1063  Double_t mMassParticle = ffP[6]*ffP[6] - (ffP[3]*ffP[3] + ffP[4]*ffP[4] + ffP[5]*ffP[5]);
1064  Double_t mMassDaughter = m[6]*m[6] - (m[3]*m[3] + m[4]*m[4] + m[5]*m[5]);
1065  if(mMassParticle > 0) mMassParticle = TMath::Sqrt(mMassParticle);
1066  if(mMassDaughter > 0) mMassDaughter = TMath::Sqrt(mMassDaughter);
1067 
1068  if( fMassHypo > -0.5)
1069  SetMassConstraint(ffP,ffC,mJ1,fMassHypo);
1070  else if((mMassParticle < SumDaughterMass) || (ffP[6]<0) )
1071  SetMassConstraint(ffP,ffC,mJ1,SumDaughterMass);
1072 
1073  if(Daughter.fMassHypo > -0.5)
1074  SetMassConstraint(m,mV,mJ2,Daughter.fMassHypo);
1075  else if((mMassDaughter < Daughter.SumDaughterMass) || (m[6] < 0) )
1076  SetMassConstraint(m,mV,mJ2,Daughter.SumDaughterMass);
1077 
1078  Double_t mDJ[7][7];
1079 
1080  for(Int_t i=0; i<7; i++) {
1081  for(Int_t j=0; j<7; j++) {
1082  mDJ[i][j] = 0;
1083  for(Int_t k=0; k<7; k++) {
1084  mDJ[i][j] += mDf[i][k]*mJ1[j][k];
1085  }
1086  }
1087  }
1088 
1089  for(Int_t i=0; i<7; ++i){
1090  for(Int_t j=0; j<7; ++j){
1091  mDf[i][j]=0;
1092  for(Int_t l=0; l<7; l++){
1093  mDf[i][j] += mJ2[i][l]*mDJ[l][j];
1094  }
1095  }
1096  }
1097 
1098  //* Add the daughter momentum to the particle momentum
1099 
1100  ffP[ 3] += m[ 3];
1101  ffP[ 4] += m[ 4];
1102  ffP[ 5] += m[ 5];
1103  ffP[ 6] += m[ 6];
1104 
1105  ffC[ 9] += mV[ 9];
1106  ffC[13] += mV[13];
1107  ffC[14] += mV[14];
1108  ffC[18] += mV[18];
1109  ffC[19] += mV[19];
1110  ffC[20] += mV[20];
1111  ffC[24] += mV[24];
1112  ffC[25] += mV[25];
1113  ffC[26] += mV[26];
1114  ffC[27] += mV[27];
1115 
1116  ffC[6 ] += mDf[3][0]; ffC[7 ] += mDf[3][1]; ffC[8 ] += mDf[3][2];
1117  ffC[10] += mDf[4][0]; ffC[11] += mDf[4][1]; ffC[12] += mDf[4][2];
1118  ffC[15] += mDf[5][0]; ffC[16] += mDf[5][1]; ffC[17] += mDf[5][2];
1119  ffC[21] += mDf[6][0]; ffC[22] += mDf[6][1]; ffC[23] += mDf[6][2];
1120 
1121  ffC[9 ] += mDf[3][3] + mDf[3][3];
1122  ffC[13] += mDf[4][3] + mDf[3][4]; ffC[14] += mDf[4][4] + mDf[4][4];
1123  ffC[18] += mDf[5][3] + mDf[3][5]; ffC[19] += mDf[5][4] + mDf[4][5]; ffC[20] += mDf[5][5] + mDf[5][5];
1124  ffC[24] += mDf[6][3] + mDf[3][6]; ffC[25] += mDf[6][4] + mDf[4][6]; ffC[26] += mDf[6][5] + mDf[5][6]; ffC[27] += mDf[6][6] + mDf[6][6];
1125 
1126  //* New estimation of the vertex position r += K*zeta
1127 
1128  for(Int_t i=0;i<7;++i)
1129  fP[i] = ffP[i];
1130 
1131  //* New covariance matrix C -= K*(mCH')'
1132 
1133  for(Int_t i=0, k=0;i<7;++i){
1134  for(Int_t j=0;j<=i;++j,++k){
1135  fC[k] = ffC[k];
1136  }
1137  }
1138  //* Calculate Chi^2
1139 
1140  fNDF += 2;
1141  fQ += Daughter.GetQ();
1142  fSFromDecay = 0;
1143  fChi2 += (mS[0]*zeta[0] + mS[1]*zeta[1] + mS[3]*zeta[2])*zeta[0]
1144  + (mS[1]*zeta[0] + mS[2]*zeta[1] + mS[4]*zeta[2])*zeta[1]
1145  + (mS[3]*zeta[0] + mS[4]*zeta[1] + mS[5]*zeta[2])*zeta[2];
1146  }
1147 }
1148 
1150 {
1151  //* Set production vertex for the particle, when the particle was not used in the vertex fit
1152 
1153  const Double_t *m = Vtx.fP, *mV = Vtx.fC;
1154 
1155  Bool_t noS = ( fC[35]<=0 ); // no decay length allowed
1156 
1157  if( noS ){
1159  fP[7] = 0;
1160  fC[28] = fC[29] = fC[30] = fC[31] = fC[32] = fC[33] = fC[34] = fC[35] = 0;
1161  } else {
1162  TransportToDS( GetDStoPoint( m ) );
1163  fP[7] = -fSFromDecay;
1164  fC[28] = fC[29] = fC[30] = fC[31] = fC[32] = fC[33] = fC[34] = 0;
1165  fC[35] = 0.1;
1166 
1167  Convert(1);
1168  }
1169 
1170  Double_t mAi[6];
1171 
1172  InvertSym3( fC, mAi );
1173 
1174  Double_t mB[5][3];
1175 
1176  mB[0][0] = fC[ 6]*mAi[0] + fC[ 7]*mAi[1] + fC[ 8]*mAi[3];
1177  mB[0][1] = fC[ 6]*mAi[1] + fC[ 7]*mAi[2] + fC[ 8]*mAi[4];
1178  mB[0][2] = fC[ 6]*mAi[3] + fC[ 7]*mAi[4] + fC[ 8]*mAi[5];
1179 
1180  mB[1][0] = fC[10]*mAi[0] + fC[11]*mAi[1] + fC[12]*mAi[3];
1181  mB[1][1] = fC[10]*mAi[1] + fC[11]*mAi[2] + fC[12]*mAi[4];
1182  mB[1][2] = fC[10]*mAi[3] + fC[11]*mAi[4] + fC[12]*mAi[5];
1183 
1184  mB[2][0] = fC[15]*mAi[0] + fC[16]*mAi[1] + fC[17]*mAi[3];
1185  mB[2][1] = fC[15]*mAi[1] + fC[16]*mAi[2] + fC[17]*mAi[4];
1186  mB[2][2] = fC[15]*mAi[3] + fC[16]*mAi[4] + fC[17]*mAi[5];
1187 
1188  mB[3][0] = fC[21]*mAi[0] + fC[22]*mAi[1] + fC[23]*mAi[3];
1189  mB[3][1] = fC[21]*mAi[1] + fC[22]*mAi[2] + fC[23]*mAi[4];
1190  mB[3][2] = fC[21]*mAi[3] + fC[22]*mAi[4] + fC[23]*mAi[5];
1191 
1192  mB[4][0] = fC[28]*mAi[0] + fC[29]*mAi[1] + fC[30]*mAi[3];
1193  mB[4][1] = fC[28]*mAi[1] + fC[29]*mAi[2] + fC[30]*mAi[4];
1194  mB[4][2] = fC[28]*mAi[3] + fC[29]*mAi[4] + fC[30]*mAi[5];
1195 
1196  Double_t z[3] = { m[0]-fP[0], m[1]-fP[1], m[2]-fP[2] };
1197 
1198  {
1199  Double_t mAVi[6] = { fC[0]-mV[0], fC[1]-mV[1], fC[2]-mV[2],
1200  fC[3]-mV[3], fC[4]-mV[4], fC[5]-mV[5] };
1201 
1202  if( !InvertSym3( mAVi, mAVi ) ){
1203 
1204  Double_t dChi2 = ( +(mAVi[0]*z[0] + mAVi[1]*z[1] + mAVi[3]*z[2])*z[0]
1205  +(mAVi[1]*z[0] + mAVi[2]*z[1] + mAVi[4]*z[2])*z[1]
1206  +(mAVi[3]*z[0] + mAVi[4]*z[1] + mAVi[5]*z[2])*z[2] );
1207 
1208  // Take Abs(dChi2) here. Negative value of 'det' or 'dChi2' shows that the particle
1209  // was not used in the production vertex fit
1210 
1211  fChi2+= TMath::Abs( dChi2 );
1212  }
1213  fNDF += 2;
1214  }
1215 
1216  fP[0] = m[0];
1217  fP[1] = m[1];
1218  fP[2] = m[2];
1219  fP[3]+= mB[0][0]*z[0] + mB[0][1]*z[1] + mB[0][2]*z[2];
1220  fP[4]+= mB[1][0]*z[0] + mB[1][1]*z[1] + mB[1][2]*z[2];
1221  fP[5]+= mB[2][0]*z[0] + mB[2][1]*z[1] + mB[2][2]*z[2];
1222  fP[6]+= mB[3][0]*z[0] + mB[3][1]*z[1] + mB[3][2]*z[2];
1223  fP[7]+= mB[4][0]*z[0] + mB[4][1]*z[1] + mB[4][2]*z[2];
1224 
1225  Double_t d0, d1, d2;
1226 
1227  fC[0] = mV[0];
1228  fC[1] = mV[1];
1229  fC[2] = mV[2];
1230  fC[3] = mV[3];
1231  fC[4] = mV[4];
1232  fC[5] = mV[5];
1233 
1234  d0= mB[0][0]*mV[0] + mB[0][1]*mV[1] + mB[0][2]*mV[3] - fC[ 6];
1235  d1= mB[0][0]*mV[1] + mB[0][1]*mV[2] + mB[0][2]*mV[4] - fC[ 7];
1236  d2= mB[0][0]*mV[3] + mB[0][1]*mV[4] + mB[0][2]*mV[5] - fC[ 8];
1237 
1238  fC[ 6]+= d0;
1239  fC[ 7]+= d1;
1240  fC[ 8]+= d2;
1241  fC[ 9]+= d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
1242 
1243  d0= mB[1][0]*mV[0] + mB[1][1]*mV[1] + mB[1][2]*mV[3] - fC[10];
1244  d1= mB[1][0]*mV[1] + mB[1][1]*mV[2] + mB[1][2]*mV[4] - fC[11];
1245  d2= mB[1][0]*mV[3] + mB[1][1]*mV[4] + mB[1][2]*mV[5] - fC[12];
1246 
1247  fC[10]+= d0;
1248  fC[11]+= d1;
1249  fC[12]+= d2;
1250  fC[13]+= d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
1251  fC[14]+= d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
1252 
1253  d0= mB[2][0]*mV[0] + mB[2][1]*mV[1] + mB[2][2]*mV[3] - fC[15];
1254  d1= mB[2][0]*mV[1] + mB[2][1]*mV[2] + mB[2][2]*mV[4] - fC[16];
1255  d2= mB[2][0]*mV[3] + mB[2][1]*mV[4] + mB[2][2]*mV[5] - fC[17];
1256 
1257  fC[15]+= d0;
1258  fC[16]+= d1;
1259  fC[17]+= d2;
1260  fC[18]+= d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
1261  fC[19]+= d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
1262  fC[20]+= d0*mB[2][0] + d1*mB[2][1] + d2*mB[2][2];
1263 
1264  d0= mB[3][0]*mV[0] + mB[3][1]*mV[1] + mB[3][2]*mV[3] - fC[21];
1265  d1= mB[3][0]*mV[1] + mB[3][1]*mV[2] + mB[3][2]*mV[4] - fC[22];
1266  d2= mB[3][0]*mV[3] + mB[3][1]*mV[4] + mB[3][2]*mV[5] - fC[23];
1267 
1268  fC[21]+= d0;
1269  fC[22]+= d1;
1270  fC[23]+= d2;
1271  fC[24]+= d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
1272  fC[25]+= d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
1273  fC[26]+= d0*mB[2][0] + d1*mB[2][1] + d2*mB[2][2];
1274  fC[27]+= d0*mB[3][0] + d1*mB[3][1] + d2*mB[3][2];
1275 
1276  d0= mB[4][0]*mV[0] + mB[4][1]*mV[1] + mB[4][2]*mV[3] - fC[28];
1277  d1= mB[4][0]*mV[1] + mB[4][1]*mV[2] + mB[4][2]*mV[4] - fC[29];
1278  d2= mB[4][0]*mV[3] + mB[4][1]*mV[4] + mB[4][2]*mV[5] - fC[30];
1279 
1280  fC[28]+= d0;
1281  fC[29]+= d1;
1282  fC[30]+= d2;
1283  fC[31]+= d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
1284  fC[32]+= d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
1285  fC[33]+= d0*mB[2][0] + d1*mB[2][1] + d2*mB[2][2];
1286  fC[34]+= d0*mB[3][0] + d1*mB[3][1] + d2*mB[3][2];
1287  fC[35]+= d0*mB[4][0] + d1*mB[4][1] + d2*mB[4][2];
1288 
1289  if( noS ){
1290  fP[7] = 0;
1291  fC[28] = fC[29] = fC[30] = fC[31] = fC[32] = fC[33] = fC[34] = fC[35] = 0;
1292  } else {
1293  TransportToDS( fP[7] );
1294  Convert(0);
1295  }
1296 
1297  fSFromDecay = 0;
1298 }
1299 
1300 void AliKFParticleBase::SetMassConstraint( Double_t *mP, Double_t *mC, Double_t mJ[7][7], Double_t mass )
1301 {
1302  //* Set nonlinear mass constraint (Mass) on the state vector mP with a covariance matrix mC.
1303 
1304  const Double_t energy2 = mP[6]*mP[6], p2 = mP[3]*mP[3]+mP[4]*mP[4]+mP[5]*mP[5], mass2 = mass*mass;
1305 
1306  const Double_t a = energy2 - p2 + 2.*mass2;
1307  const Double_t b = -2.*(energy2 + p2);
1308  const Double_t c = energy2 - p2 - mass2;
1309 
1310  Double_t lambda = 0;
1311  if(TMath::Abs(b) > 1.e-10) lambda = -c / b;
1312 
1313  Double_t d = 4.*energy2*p2 - mass2*(energy2-p2-2.*mass2);
1314  if(d>=0 && TMath::Abs(a) > 1.e-10) lambda = (energy2 + p2 - sqrt(d))/a;
1315 
1316  if(mP[6] < 0) //If energy < 0 we need a lambda < 0
1317  lambda = -1000000.; //Empirical, a better solution should be found
1318 
1319  Int_t iIter=0;
1320  for(iIter=0; iIter<100; iIter++)
1321  {
1322  Double_t lambda2 = lambda*lambda;
1323  Double_t lambda4 = lambda2*lambda2;
1324 
1325  Double_t lambda0 = lambda;
1326 
1327  Double_t f = -mass2 * lambda4 + a*lambda2 + b*lambda + c;
1328  Double_t df = -4.*mass2 * lambda2*lambda + 2.*a*lambda + b;
1329  if(TMath::Abs(df) < 1.e-10) break;
1330  lambda -= f/df;
1331  if(TMath::Abs(lambda0 - lambda) < 1.e-8) break;
1332  }
1333 
1334  const Double_t lpi = 1./(1. + lambda);
1335  const Double_t lmi = 1./(1. - lambda);
1336  const Double_t lp2i = lpi*lpi;
1337  const Double_t lm2i = lmi*lmi;
1338 
1339  Double_t lambda2 = lambda*lambda;
1340 
1341  Double_t dfl = -4.*mass2 * lambda2*lambda + 2.*a*lambda + b;
1342  Double_t dfx[7] = {0};//,0,0,0};
1343  dfx[0] = -2.*(1. + lambda)*(1. + lambda)*mP[3];
1344  dfx[1] = -2.*(1. + lambda)*(1. + lambda)*mP[4];
1345  dfx[2] = -2.*(1. + lambda)*(1. + lambda)*mP[5];
1346  dfx[3] = 2.*(1. - lambda)*(1. - lambda)*mP[6];
1347  Double_t dlx[4] = {1,1,1,1};
1348  if(TMath::Abs(dfl) > 1.e-10 )
1349  {
1350  for(int i=0; i<4; i++)
1351  dlx[i] = -dfx[i] / dfl;
1352  }
1353 
1354  Double_t dxx[4] = {mP[3]*lm2i, mP[4]*lm2i, mP[5]*lm2i, -mP[6]*lp2i};
1355 
1356  for(Int_t i=0; i<7; i++)
1357  for(Int_t j=0; j<7; j++)
1358  mJ[i][j]=0;
1359  mJ[0][0] = 1.;
1360  mJ[1][1] = 1.;
1361  mJ[2][2] = 1.;
1362 
1363  for(Int_t i=3; i<7; i++)
1364  for(Int_t j=3; j<7; j++)
1365  mJ[i][j] = dlx[j-3]*dxx[i-3];
1366 
1367  for(Int_t i=3; i<6; i++)
1368  mJ[i][i] += lmi;
1369  mJ[6][6] += lpi;
1370 
1371  Double_t mCJ[7][7];
1372 
1373  for(Int_t i=0; i<7; i++) {
1374  for(Int_t j=0; j<7; j++) {
1375  mCJ[i][j] = 0;
1376  for(Int_t k=0; k<7; k++) {
1377  mCJ[i][j] += mC[IJ(i,k)]*mJ[j][k];
1378  }
1379  }
1380  }
1381 
1382  for(Int_t i=0; i<7; ++i){
1383  for(Int_t j=0; j<=i; ++j){
1384  mC[IJ(i,j)]=0;
1385  for(Int_t l=0; l<7; l++){
1386  mC[IJ(i,j)] += mJ[i][l]*mCJ[l][j];
1387  }
1388  }
1389  }
1390 
1391  mP[3] *= lmi;
1392  mP[4] *= lmi;
1393  mP[5] *= lmi;
1394  mP[6] *= lpi;
1395 }
1396 
1398 {
1399  //* Set nonlinear mass constraint (mass)
1400 
1401  Double_t mJ[7][7];
1402  SetMassConstraint( fP, fC, mJ, mass );
1403  fMassHypo = mass;
1404  SumDaughterMass = mass;
1405 }
1406 
1407 void AliKFParticleBase::SetMassConstraint( Double_t Mass, Double_t SigmaMass )
1408 {
1409  //* Set hard( SigmaMass=0 ) or soft (SigmaMass>0) mass constraint
1410 
1411  fMassHypo = Mass;
1412  SumDaughterMass = Mass;
1413 
1414  Double_t m2 = Mass*Mass; // measurement, weighted by Mass
1415  Double_t s2 = m2*SigmaMass*SigmaMass; // sigma^2
1416 
1417  Double_t p2 = fP[3]*fP[3] + fP[4]*fP[4] + fP[5]*fP[5];
1418  Double_t e0 = TMath::Sqrt(m2+p2);
1419 
1420  Double_t mH[8];
1421  mH[0] = mH[1] = mH[2] = 0.;
1422  mH[3] = -2*fP[3];
1423  mH[4] = -2*fP[4];
1424  mH[5] = -2*fP[5];
1425  mH[6] = 2*fP[6];//e0;
1426  mH[7] = 0;
1427 
1428  Double_t zeta = e0*e0 - e0*fP[6];
1429  zeta = m2 - (fP[6]*fP[6]-p2);
1430 
1431  Double_t mCHt[8], s2_est=0;
1432  for( Int_t i=0; i<8; ++i ){
1433  mCHt[i] = 0.0;
1434  for (Int_t j=0;j<8;++j) mCHt[i] += Cij(i,j)*mH[j];
1435  s2_est += mH[i]*mCHt[i];
1436  }
1437 
1438  if( s2_est<1.e-20 ) return; // calculated mass error is already 0,
1439  // the particle can not be constrained on mass
1440 
1441  Double_t w2 = 1./( s2 + s2_est );
1442  fChi2 += zeta*zeta*w2;
1443  fNDF += 1;
1444  for( Int_t i=0, ii=0; i<8; ++i ){
1445  Double_t ki = mCHt[i]*w2;
1446  fP[i]+= ki*zeta;
1447  for(Int_t j=0;j<=i;++j) fC[ii++] -= ki*mCHt[j];
1448  }
1449 }
1450 
1451 
1453 {
1454  //* Set no decay length for resonances
1455 
1457 
1458  Double_t h[8];
1459  h[0] = h[1] = h[2] = h[3] = h[4] = h[5] = h[6] = 0;
1460  h[7] = 1;
1461 
1462  Double_t zeta = 0 - fP[7];
1463  for(Int_t i=0;i<8;++i) zeta -= h[i]*(fP[i]-fP[i]);
1464 
1465  Double_t s = fC[35];
1466  if( s>1.e-20 ){
1467  s = 1./s;
1468  fChi2 += zeta*zeta*s;
1469  fNDF += 1;
1470  for( Int_t i=0, ii=0; i<7; ++i ){
1471  Double_t ki = fC[28+i]*s;
1472  fP[i]+= ki*zeta;
1473  for(Int_t j=0;j<=i;++j) fC[ii++] -= ki*fC[28+j];
1474  }
1475  }
1476  fP[7] = 0;
1477  fC[28] = fC[29] = fC[30] = fC[31] = fC[32] = fC[33] = fC[34] = fC[35] = 0;
1478 }
1479 
1480 
1481 void AliKFParticleBase::Construct( const AliKFParticleBase* vDaughters[], Int_t NDaughters,
1482  const AliKFParticleBase *Parent, Double_t Mass, Bool_t IsConstrained )
1483 {
1484  //* Full reconstruction in one go
1485 
1486  Int_t maxIter = 1;
1487  bool wasLinearized = fIsLinearized;
1488  if( !fIsLinearized || IsConstrained ){
1489  //fVtxGuess[0] = fVtxGuess[1] = fVtxGuess[2] = 0; //!!!!
1490  fVtxGuess[0] = GetX();
1491  fVtxGuess[1] = GetY();
1492  fVtxGuess[2] = GetZ();
1493  fIsLinearized = 1;
1494  maxIter = 3;
1495  }
1496 
1497  Double_t constraintC[6];
1498 
1499  if( IsConstrained ){
1500  for(Int_t i=0;i<6;++i) constraintC[i]=fC[i];
1501  } else {
1502  for(Int_t i=0;i<6;++i) constraintC[i]=0.;
1503  constraintC[0] = constraintC[2] = constraintC[5] = 100.;
1504  }
1505 
1506 
1507  for( Int_t iter=0; iter<maxIter; iter++ ){
1508  fAtProductionVertex = 0;
1509  fSFromDecay = 0;
1510  fP[0] = fVtxGuess[0];
1511  fP[1] = fVtxGuess[1];
1512  fP[2] = fVtxGuess[2];
1513  fP[3] = 0;
1514  fP[4] = 0;
1515  fP[5] = 0;
1516  fP[6] = 0;
1517  fP[7] = 0;
1518  SumDaughterMass = 0;
1519 
1520  for(Int_t i=0;i<6; ++i) fC[i]=constraintC[i];
1521  for(Int_t i=6;i<36;++i) fC[i]=0.;
1522  fC[35] = 1.;
1523 
1524  fNDF = IsConstrained ?0 :-3;
1525  fChi2 = 0.;
1526  fQ = 0;
1527 
1528  for( Int_t itr =0; itr<NDaughters; itr++ ){
1529  AddDaughter( *vDaughters[itr] );
1530  }
1531  if( iter<maxIter-1){
1532  for( Int_t i=0; i<3; i++ ) fVtxGuess[i] = fP[i];
1533  }
1534  }
1535  fIsLinearized = wasLinearized;
1536 
1537  if( Mass>=0 ) SetMassConstraint( Mass );
1538  if( Parent ) SetProductionVertex( *Parent );
1539 }
1540 
1541 
1542 void AliKFParticleBase::Convert( bool ToProduction )
1543 {
1544  //* Tricky function - convert the particle error along its trajectory to
1545  //* the value which corresponds to its production/decay vertex
1546  //* It is done by combination of the error of decay length with the position errors
1547 
1548  Double_t fld[3];
1549  {
1550  GetFieldValue( fP, fld );
1551  const Double_t kCLight = fQ*0.000299792458;
1552  fld[0]*=kCLight; fld[1]*=kCLight; fld[2]*=kCLight;
1553  }
1554 
1555  Double_t h[6];
1556 
1557  h[0] = fP[3];
1558  h[1] = fP[4];
1559  h[2] = fP[5];
1560  if( ToProduction ){ h[0]=-h[0]; h[1]=-h[1]; h[2]=-h[2]; }
1561  h[3] = h[1]*fld[2]-h[2]*fld[1];
1562  h[4] = h[2]*fld[0]-h[0]*fld[2];
1563  h[5] = h[0]*fld[1]-h[1]*fld[0];
1564 
1565  Double_t c;
1566 
1567  c = fC[28]+h[0]*fC[35];
1568  fC[ 0]+= h[0]*(c+fC[28]);
1569  fC[28] = c;
1570 
1571  fC[ 1]+= h[1]*fC[28] + h[0]*fC[29];
1572  c = fC[29]+h[1]*fC[35];
1573  fC[ 2]+= h[1]*(c+fC[29]);
1574  fC[29] = c;
1575 
1576  fC[ 3]+= h[2]*fC[28] + h[0]*fC[30];
1577  fC[ 4]+= h[2]*fC[29] + h[1]*fC[30];
1578  c = fC[30]+h[2]*fC[35];
1579  fC[ 5]+= h[2]*(c+fC[30]);
1580  fC[30] = c;
1581 
1582  fC[ 6]+= h[3]*fC[28] + h[0]*fC[31];
1583  fC[ 7]+= h[3]*fC[29] + h[1]*fC[31];
1584  fC[ 8]+= h[3]*fC[30] + h[2]*fC[31];
1585  c = fC[31]+h[3]*fC[35];
1586  fC[ 9]+= h[3]*(c+fC[31]);
1587  fC[31] = c;
1588 
1589  fC[10]+= h[4]*fC[28] + h[0]*fC[32];
1590  fC[11]+= h[4]*fC[29] + h[1]*fC[32];
1591  fC[12]+= h[4]*fC[30] + h[2]*fC[32];
1592  fC[13]+= h[4]*fC[31] + h[3]*fC[32];
1593  c = fC[32]+h[4]*fC[35];
1594  fC[14]+= h[4]*(c+fC[32]);
1595  fC[32] = c;
1596 
1597  fC[15]+= h[5]*fC[28] + h[0]*fC[33];
1598  fC[16]+= h[5]*fC[29] + h[1]*fC[33];
1599  fC[17]+= h[5]*fC[30] + h[2]*fC[33];
1600  fC[18]+= h[5]*fC[31] + h[3]*fC[33];
1601  fC[19]+= h[5]*fC[32] + h[4]*fC[33];
1602  c = fC[33]+h[5]*fC[35];
1603  fC[20]+= h[5]*(c+fC[33]);
1604  fC[33] = c;
1605 
1606  fC[21]+= h[0]*fC[34];
1607  fC[22]+= h[1]*fC[34];
1608  fC[23]+= h[2]*fC[34];
1609  fC[24]+= h[3]*fC[34];
1610  fC[25]+= h[4]*fC[34];
1611  fC[26]+= h[5]*fC[34];
1612 }
1613 
1614 
1616 {
1617  //* Transport the particle to its decay vertex
1618 
1619  if( fSFromDecay != 0 ) TransportToDS( -fSFromDecay );
1620  if( fAtProductionVertex ) Convert(0);
1621  fAtProductionVertex = 0;
1622 }
1623 
1625 {
1626  //* Transport the particle to its production vertex
1627 
1628  if( fSFromDecay != -fP[7] ) TransportToDS( -fSFromDecay-fP[7] );
1629  if( !fAtProductionVertex ) Convert( 1 );
1630  fAtProductionVertex = 1;
1631 }
1632 
1633 
1635 {
1636  //* Transport the particle on dS parameter (SignedPath/Momentum)
1637 
1638  Transport( dS, fP, fC );
1639  fSFromDecay+= dS;
1640 }
1641 
1642 
1643 Double_t AliKFParticleBase::GetDStoPointLine( const Double_t xyz[] ) const
1644 {
1645  //* Get dS to a certain space point without field
1646 
1647  Double_t p2 = fP[3]*fP[3] + fP[4]*fP[4] + fP[5]*fP[5];
1648  if( p2<1.e-4 ) p2 = 1;
1649  return ( fP[3]*(xyz[0]-fP[0]) + fP[4]*(xyz[1]-fP[1]) + fP[5]*(xyz[2]-fP[2]) )/p2;
1650 }
1651 
1652 
1653 Double_t AliKFParticleBase::GetDStoPointBz( Double_t B, const Double_t xyz[] )
1654  const
1655 {
1656 
1657  //* Get dS to a certain space point for Bz field
1658  const Double_t kCLight = 0.000299792458;
1659  Double_t bq = B*fQ*kCLight;
1660  Double_t pt2 = fP[3]*fP[3] + fP[4]*fP[4];
1661  if( pt2<1.e-4 ) return 0;
1662  Double_t dx = xyz[0] - fP[0];
1663  Double_t dy = xyz[1] - fP[1];
1664  Double_t a = dx*fP[3]+dy*fP[4];
1665  Double_t dS;
1666 
1667  if( TMath::Abs(bq)<1.e-8 ) dS = a/pt2;
1668  else dS = TMath::ATan2( bq*a, pt2 + bq*(dy*fP[3] -dx*fP[4]) )/bq;
1669 
1670  if(0){
1671 
1672  Double_t px = fP[3];
1673  Double_t py = fP[4];
1674  Double_t pz = fP[5];
1675  Double_t ss[2], g[2][5];
1676 
1677  ss[0] = dS;
1678  ss[1] = -dS;
1679  for( Int_t i=0; i<2; i++){
1680  Double_t bs = bq*ss[i];
1681  Double_t c = TMath::Cos(bs), s = TMath::Sin(bs);
1682  Double_t cB,sB;
1683  if( TMath::Abs(bq)>1.e-8){
1684  cB= (1-c)/bq;
1685  sB= s/bq;
1686  }else{
1687  const Double_t kOvSqr6 = 1./TMath::Sqrt(6.);
1688  sB = (1.-bs*kOvSqr6)*(1.+bs*kOvSqr6)*ss[i];
1689  cB = .5*sB*bs;
1690  }
1691  g[i][0] = fP[0] + sB*px + cB*py;
1692  g[i][1] = fP[1] - cB*px + sB*py;
1693  g[i][2] = fP[2] + ss[i]*pz;
1694  g[i][3] = + c*px + s*py;
1695  g[i][4] = - s*px + c*py;
1696  }
1697 
1698  Int_t i=0;
1699 
1700  Double_t dMin = 1.e10;
1701  for( Int_t j=0; j<2; j++){
1702  Double_t xx = g[j][0]-xyz[0];
1703  Double_t yy = g[j][1]-xyz[1];
1704  Double_t zz = g[j][2]-xyz[2];
1705  Double_t d = xx*xx + yy*yy + zz*zz;
1706  if( d<dMin ){
1707  dMin = d;
1708  i = j;
1709  }
1710  }
1711 
1712  dS = ss[i];
1713 
1714  Double_t x= g[i][0], y= g[i][1], z= g[i][2], ppx= g[i][3], ppy= g[i][4];
1715  Double_t ddx = x-xyz[0];
1716  Double_t ddy = y-xyz[1];
1717  Double_t ddz = z-xyz[2];
1718  Double_t c = ddx*ppx + ddy*ppy + ddz*pz ;
1719  Double_t pp2 = ppx*ppx + ppy*ppy + pz*pz;
1720  if( TMath::Abs(pp2)>1.e-8 ){
1721  dS+=c/pp2;
1722  }
1723  }
1724  return dS;
1725 }
1726 
1727 
1729  Double_t &DS, Double_t &DS1 )
1730  const
1731 {
1732  //* Get dS to another particle for Bz field
1733  Double_t px = fP[3];
1734  Double_t py = fP[4];
1735  Double_t pz = fP[5];
1736 
1737  Double_t px1 = p.fP[3];
1738  Double_t py1 = p.fP[4];
1739  Double_t pz1 = p.fP[5];
1740 
1741  const Double_t kCLight = 0.000299792458;
1742 
1743  Double_t bq = B*fQ*kCLight;
1744  Double_t bq1 = B*p.fQ*kCLight;
1745  Double_t s=0, ds=0, s1=0, ds1=0;
1746 
1747  if( TMath::Abs(bq)>1.e-8 || TMath::Abs(bq1)>1.e-8 ){
1748 
1749  Double_t dx = (p.fP[0] - fP[0]);
1750  Double_t dy = (p.fP[1] - fP[1]);
1751  Double_t d2 = (dx*dx+dy*dy);
1752 
1753  Double_t p2 = (px *px + py *py);
1754  Double_t p21 = (px1*px1 + py1*py1);
1755 
1756  if( TMath::Abs(p2) < 1.e-8 || TMath::Abs(p21) < 1.e-8 )
1757  {
1758  DS=0.;
1759  DS1=0.;
1760  return;
1761  }
1762 
1763  Double_t a = (px*py1 - py*px1);
1764  Double_t b = (px*px1 + py*py1);
1765 
1766  Double_t ldx = bq*bq1*dx - bq1*py + bq*py1 ;
1767  Double_t ldy = bq*bq1*dy + bq1*px - bq*px1 ;
1768  Double_t l2 = ldx*ldx + ldy*ldy;
1769 
1770  Double_t cS = bq1*p2 + bq*bq1*(dy* px - dx* py) - bq*b;
1771  Double_t cS1= bq*p21 - bq*bq1*(dy*px1 - dx*py1) - bq1*b;
1772 
1773  Double_t ca = bq*bq*bq1*d2 +2*( cS + bq*bq*(py1*dx-px1*dy)) ;
1774  Double_t ca1 = bq*bq1*bq1*d2 +2*( cS1 - bq1*bq1*(py*dx-px*dy)) ;
1775 
1776  Double_t sa = 4*l2*p2 - ca*ca;
1777  Double_t sa1 = 4*l2*p21 - ca1*ca1;
1778 
1779  if(sa<0) sa=0;
1780  if(sa1<0)sa1=0;
1781 
1782  if( TMath::Abs(bq)>1.e-8){
1783  s = TMath::ATan2( bq*( bq1*(dx*px +dy*py) + a ) , cS )/bq;
1784  ds = TMath::ATan2(TMath::Sqrt(sa),ca)/bq;
1785  } else {
1786  s = ( (dx*px + dy*py) + (py*px1-px*py1)/bq1)/p2;
1787  ds = s*s - (d2-2*(px1*dy-py1*dx)/bq1)/p2;
1788  if( ds<0 ) ds = 0;
1789  ds = TMath::Sqrt(ds);
1790  }
1791 
1792  if( TMath::Abs(bq1)>1.e-8){
1793  s1 = TMath::ATan2( -bq1*( bq*(dx*px1+dy*py1) + a), cS1 )/bq1;
1794  ds1 = TMath::ATan2(TMath::Sqrt(sa1),ca1)/bq1;
1795  } else {
1796  s1 = (-(dx*px1 + dy*py1) + (py*px1-px*py1)/bq)/p21;
1797  ds1 = s1*s1 - (d2+2*(px*dy-py*dx)/bq)/p21;
1798  if( ds1<0 ) ds1 = 0;
1799  ds1 = TMath::Sqrt(ds1);
1800  }
1801  }
1802 
1803  Double_t ss[2], ss1[2], g[2][5],g1[2][5];
1804 
1805  ss[0] = s + ds;
1806  ss[1] = s - ds;
1807  ss1[0] = s1 + ds1;
1808  ss1[1] = s1 - ds1;
1809  for( Int_t i=0; i<2; i++){
1810  Double_t bs = bq*ss[i];
1811  Double_t c = TMath::Cos(bs), sss = TMath::Sin(bs);
1812  Double_t cB,sB;
1813  if( TMath::Abs(bq)>1.e-8){
1814  cB= (1-c)/bq;
1815  sB= sss/bq;
1816  }else{
1817  const Double_t kOvSqr6 = 1./TMath::Sqrt(6.);
1818  sB = (1.-bs*kOvSqr6)*(1.+bs*kOvSqr6)*ss[i];
1819  cB = .5*sB*bs;
1820  }
1821  g[i][0] = fP[0] + sB*px + cB*py;
1822  g[i][1] = fP[1] - cB*px + sB*py;
1823  g[i][2] = fP[2] + ss[i]*pz;
1824  g[i][3] = + c*px + sss*py;
1825  g[i][4] = - sss*px + c*py;
1826 
1827  bs = bq1*ss1[i];
1828  c = TMath::Cos(bs); sss = TMath::Sin(bs);
1829  if( TMath::Abs(bq1)>1.e-8){
1830  cB= (1-c)/bq1;
1831  sB= sss/bq1;
1832  }else{
1833  const Double_t kOvSqr6 = 1./TMath::Sqrt(6.);
1834  sB = (1.-bs*kOvSqr6)*(1.+bs*kOvSqr6)*ss1[i];
1835  cB = .5*sB*bs;
1836  }
1837 
1838  g1[i][0] = p.fP[0] + sB*px1 + cB*py1;
1839  g1[i][1] = p.fP[1] - cB*px1 + sB*py1;
1840  g1[i][2] = p.fP[2] + ss[i]*pz1;
1841  g1[i][3] = + c*px1 + sss*py1;
1842  g1[i][4] = - sss*px1 + c*py1;
1843  }
1844 
1845  Int_t i=0, i1=0;
1846 
1847  Double_t dMin = 1.e10;
1848  for( Int_t j=0; j<2; j++){
1849  for( Int_t j1=0; j1<2; j1++){
1850  Double_t xx = g[j][0]-g1[j1][0];
1851  Double_t yy = g[j][1]-g1[j1][1];
1852  Double_t zz = g[j][2]-g1[j1][2];
1853  Double_t d = xx*xx + yy*yy + zz*zz;
1854  if( d<dMin ){
1855  dMin = d;
1856  i = j;
1857  i1 = j1;
1858  }
1859  }
1860  }
1861 
1862  DS = ss[i];
1863  DS1 = ss1[i1];
1864  if(0){
1865  Double_t x= g[i][0], y= g[i][1], z= g[i][2], ppx= g[i][3], ppy= g[i][4];
1866  Double_t x1=g1[i1][0], y1= g1[i1][1], z1= g1[i1][2], ppx1= g1[i1][3], ppy1= g1[i1][4];
1867  Double_t dx = x1-x;
1868  Double_t dy = y1-y;
1869  Double_t dz = z1-z;
1870  Double_t a = ppx*ppx1 + ppy*ppy1 + pz*pz1;
1871  Double_t b = dx*ppx1 + dy*ppy1 + dz*pz1;
1872  Double_t c = dx*ppx + dy*ppy + dz*pz ;
1873  Double_t pp2 = ppx*ppx + ppy*ppy + pz*pz;
1874  Double_t pp21= ppx1*ppx1 + ppy1*ppy1 + pz1*pz1;
1875  Double_t det = pp2*pp21 - a*a;
1876  if( TMath::Abs(det)>1.e-8 ){
1877  DS+=(a*b-pp21*c)/det;
1878  DS1+=(a*c-pp2*b)/det;
1879  }
1880  }
1881 }
1882 
1883 
1884 
1886  Double_t P[], Double_t C[] ) const
1887 {
1888  //* Transport the particle on dS, output to P[],C[], for CBM field
1889 
1890  if( fQ==0 ){
1891  TransportLine( dS, P, C );
1892  return;
1893  }
1894 
1895  const Double_t kCLight = 0.000299792458;
1896 
1897  Double_t c = fQ*kCLight;
1898 
1899  // construct coefficients
1900 
1901  Double_t
1902  px = fP[3],
1903  py = fP[4],
1904  pz = fP[5];
1905 
1906  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;
1907 
1908  { // get field integrals
1909 
1910  Double_t fld[3][3];
1911  Double_t p0[3], p1[3], p2[3];
1912 
1913  // line track approximation
1914 
1915  p0[0] = fP[0];
1916  p0[1] = fP[1];
1917  p0[2] = fP[2];
1918 
1919  p2[0] = fP[0] + px*dS;
1920  p2[1] = fP[1] + py*dS;
1921  p2[2] = fP[2] + pz*dS;
1922 
1923  p1[0] = 0.5*(p0[0]+p2[0]);
1924  p1[1] = 0.5*(p0[1]+p2[1]);
1925  p1[2] = 0.5*(p0[2]+p2[2]);
1926 
1927  // first order track approximation
1928  {
1929  GetFieldValue( p0, fld[0] );
1930  GetFieldValue( p1, fld[1] );
1931  GetFieldValue( p2, fld[2] );
1932 
1933  Double_t ssy1 = ( 7*fld[0][1] + 6*fld[1][1]-fld[2][1] )*c*dS*dS/96.;
1934  Double_t ssy2 = ( fld[0][1] + 2*fld[1][1] )*c*dS*dS/6.;
1935 
1936  p1[0] -= ssy1*pz;
1937  p1[2] += ssy1*px;
1938  p2[0] -= ssy2*pz;
1939  p2[2] += ssy2*px;
1940  }
1941 
1942  GetFieldValue( p0, fld[0] );
1943  GetFieldValue( p1, fld[1] );
1944  GetFieldValue( p2, fld[2] );
1945 
1946  sx = c*( fld[0][0] + 4*fld[1][0] + fld[2][0] )*dS/6.;
1947  sy = c*( fld[0][1] + 4*fld[1][1] + fld[2][1] )*dS/6.;
1948  sz = c*( fld[0][2] + 4*fld[1][2] + fld[2][2] )*dS/6.;
1949 
1950  ssx = c*( fld[0][0] + 2*fld[1][0])*dS*dS/6.;
1951  ssy = c*( fld[0][1] + 2*fld[1][1])*dS*dS/6.;
1952  ssz = c*( fld[0][2] + 2*fld[1][2])*dS*dS/6.;
1953 
1954  Double_t c2[3][3] = { { 5, -4, -1},{ 44, 80, -4},{ 11, 44, 5} }; // /=360.
1955  Double_t cc2[3][3] = { { 38, 8, -4},{ 148, 208, -20},{ 3, 36, 3} }; // /=2520.
1956  for(Int_t n=0; n<3; n++)
1957  for(Int_t m=0; m<3; m++)
1958  {
1959  syz += c2[n][m]*fld[n][1]*fld[m][2];
1960  ssyz += cc2[n][m]*fld[n][1]*fld[m][2];
1961  }
1962 
1963  syz *= c*c*dS*dS/360.;
1964  ssyz *= c*c*dS*dS*dS/2520.;
1965 
1966  syy = c*( fld[0][1] + 4*fld[1][1] + fld[2][1] )*dS;
1967  syyy = syy*syy*syy / 1296;
1968  syy = syy*syy/72;
1969 
1970  ssyy = ( fld[0][1]*( 38*fld[0][1] + 156*fld[1][1] - fld[2][1] )+
1971  fld[1][1]*( 208*fld[1][1] +16*fld[2][1] )+
1972  fld[2][1]*( 3*fld[2][1] )
1973  )*dS*dS*dS*c*c/2520.;
1974  ssyyy =
1975  (
1976  fld[0][1]*( fld[0][1]*( 85*fld[0][1] + 526*fld[1][1] - 7*fld[2][1] )+
1977  fld[1][1]*( 1376*fld[1][1] +84*fld[2][1] )+
1978  fld[2][1]*( 19*fld[2][1] ) )+
1979  fld[1][1]*( fld[1][1]*( 1376*fld[1][1] +256*fld[2][1] )+
1980  fld[2][1]*( 62*fld[2][1] ) )+
1981  fld[2][1]*fld[2][1] *( 3*fld[2][1] )
1982  )*dS*dS*dS*dS*c*c*c/90720.;
1983 
1984  }
1985 
1986  Double_t mJ[8][8];
1987  for( Int_t i=0; i<8; i++ ) for( Int_t j=0; j<8; j++) mJ[i][j]=0;
1988 
1989  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;
1990  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;
1991  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;
1992 
1993  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;
1994  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;
1995  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;
1996  mJ[6][6] = mJ[7][7] = 1;
1997 
1998  P[0] = fP[0] + mJ[0][3]*px + mJ[0][4]*py + mJ[0][5]*pz;
1999  P[1] = fP[1] + mJ[1][3]*px + mJ[1][4]*py + mJ[1][5]*pz;
2000  P[2] = fP[2] + mJ[2][3]*px + mJ[2][4]*py + mJ[2][5]*pz;
2001  P[3] = mJ[3][3]*px + mJ[3][4]*py + mJ[3][5]*pz;
2002  P[4] = mJ[4][3]*px + mJ[4][4]*py + mJ[4][5]*pz;
2003  P[5] = mJ[5][3]*px + mJ[5][4]*py + mJ[5][5]*pz;
2004  P[6] = fP[6];
2005  P[7] = fP[7];
2006 
2007  MultQSQt( mJ[0], fC, C);
2008 
2009 }
2010 
2011 
2012 void AliKFParticleBase::TransportBz( Double_t b, Double_t t,
2013  Double_t p[], Double_t e[] ) const
2014 {
2015  //* Transport the particle on dS, output to P[],C[], for Bz field
2016 
2017  const Double_t kCLight = 0.000299792458;
2018  b = b*fQ*kCLight;
2019  Double_t bs= b*t;
2020  Double_t s = TMath::Sin(bs), c = TMath::Cos(bs);
2021  Double_t sB, cB;
2022  if( TMath::Abs(bs)>1.e-10){
2023  sB= s/b;
2024  cB= (1-c)/b;
2025  }else{
2026  const Double_t kOvSqr6 = 1./TMath::Sqrt(6.);
2027  sB = (1.-bs*kOvSqr6)*(1.+bs*kOvSqr6)*t;
2028  cB = .5*sB*bs;
2029  }
2030 
2031  Double_t px = fP[3];
2032  Double_t py = fP[4];
2033  Double_t pz = fP[5];
2034 
2035  p[0] = fP[0] + sB*px + cB*py;
2036  p[1] = fP[1] - cB*px + sB*py;
2037  p[2] = fP[2] + t*pz;
2038  p[3] = c*px + s*py;
2039  p[4] = -s*px + c*py;
2040  p[5] = fP[5];
2041  p[6] = fP[6];
2042  p[7] = fP[7];
2043 
2044  /*
2045  Double_t mJ[8][8] = { {1,0,0, sB, cB, 0, 0, 0 },
2046  {0,1,0, -cB, sB, 0, 0, 0 },
2047  {0,0,1, 0, 0, t, 0, 0 },
2048  {0,0,0, c, s, 0, 0, 0 },
2049  {0,0,0, -s, c, 0, 0, 0 },
2050  {0,0,0, 0, 0, 1, 0, 0 },
2051  {0,0,0, 0, 0, 0, 1, 0 },
2052  {0,0,0, 0, 0, 0, 0, 1 } };
2053  Double_t mA[8][8];
2054  for( Int_t k=0,i=0; i<8; i++)
2055  for( Int_t j=0; j<=i; j++, k++ ) mA[i][j] = mA[j][i] = fC[k];
2056 
2057  Double_t mJC[8][8];
2058  for( Int_t i=0; i<8; i++ )
2059  for( Int_t j=0; j<8; j++ ){
2060  mJC[i][j]=0;
2061  for( Int_t k=0; k<8; k++ ) mJC[i][j]+=mJ[i][k]*mA[k][j];
2062  }
2063 
2064  for( Int_t k=0,i=0; i<8; i++)
2065  for( Int_t j=0; j<=i; j++, k++ ){
2066  e[k] = 0;
2067  for( Int_t l=0; l<8; l++ ) e[k]+=mJC[i][l]*mJ[j][l];
2068  }
2069 
2070  return;
2071  */
2072 
2073  Double_t
2074  c6=fC[6], c7=fC[7], c8=fC[8], c17=fC[17], c18=fC[18],
2075  c24 = fC[24], c31 = fC[31];
2076 
2077  Double_t
2078  cBC13 = cB*fC[13],
2079  mJC13 = c7 - cB*fC[9] + sB*fC[13],
2080  mJC14 = fC[11] - cBC13 + sB*fC[14],
2081  mJC23 = c8 + t*c18,
2082  mJC24 = fC[12] + t*fC[19],
2083  mJC33 = c*fC[9] + s*fC[13],
2084  mJC34 = c*fC[13] + s*fC[14],
2085  mJC43 = -s*fC[9] + c*fC[13],
2086  mJC44 = -s*fC[13] + c*fC[14];
2087 
2088 
2089  e[0]= fC[0] + 2*(sB*c6 + cB*fC[10]) + (sB*fC[9] + 2*cBC13)*sB + cB*cB*fC[14];
2090  e[1]= fC[1] - cB*c6 + sB*fC[10] + mJC13*sB + mJC14*cB;
2091  e[2]= fC[2] - cB*c7 + sB*fC[11] - mJC13*cB + mJC14*sB;
2092  e[3]= fC[3] + t*fC[15] + mJC23*sB + mJC24*cB;
2093  e[4]= fC[4] + t*fC[16] - mJC23*cB + mJC24*sB;
2094 
2095  e[15]= fC[15] + c18*sB + fC[19]*cB;
2096  e[16]= fC[16] - c18*cB + fC[19]*sB;
2097  e[17]= c17 + fC[20]*t;
2098  e[18]= c18*c + fC[19]*s;
2099  e[19]= -c18*s + fC[19]*c;
2100 
2101  e[5]= fC[5] + (c17 + e[17] )*t;
2102 
2103  e[6]= c*c6 + s*fC[10] + mJC33*sB + mJC34*cB;
2104  e[7]= c*c7 + s*fC[11] - mJC33*cB + mJC34*sB;
2105  e[8]= c*c8 + s*fC[12] + e[18]*t;
2106  e[9]= mJC33*c + mJC34*s;
2107  e[10]= -s*c6 + c*fC[10] + mJC43*sB + mJC44*cB;
2108 
2109 
2110  e[11]= -s*c7 + c*fC[11] - mJC43*cB + mJC44*sB;
2111  e[12]= -s*c8 + c*fC[12] + e[19]*t;
2112  e[13]= mJC43*c + mJC44*s;
2113  e[14]= -mJC43*s + mJC44*c;
2114  e[20]= fC[20];
2115  e[21]= fC[21] + fC[25]*cB + c24*sB;
2116  e[22]= fC[22] - c24*cB + fC[25]*sB;
2117  e[23]= fC[23] + fC[26]*t;
2118  e[24]= c*c24 + s*fC[25];
2119  e[25]= c*fC[25] - c24*s;
2120  e[26]= fC[26];
2121  e[27]= fC[27];
2122  e[28]= fC[28] + fC[32]*cB + c31*sB;
2123  e[29]= fC[29] - c31*cB + fC[32]*sB;
2124  e[30]= fC[30] + fC[33]*t;
2125  e[31]= c*c31 + s*fC[32];
2126  e[32]= c*fC[32] - s*c31;
2127  e[33]= fC[33];
2128  e[34]= fC[34];
2129  e[35]= fC[35];
2130 }
2131 
2132 
2134 {
2135  //* Calculate distance from vertex [cm]
2136 
2137  return GetDistanceFromVertex( Vtx.fP );
2138 }
2139 
2140 Double_t AliKFParticleBase::GetDistanceFromVertex( const Double_t vtx[] ) const
2141 {
2142  //* Calculate distance from vertex [cm]
2143 
2144  Double_t mP[8], mC[36];
2145  Transport( GetDStoPoint(vtx), mP, mC );
2146  Double_t d[3]={ vtx[0]-mP[0], vtx[1]-mP[1], vtx[2]-mP[2]};
2147  return TMath::Sqrt( d[0]*d[0]+d[1]*d[1]+d[2]*d[2] );
2148 }
2149 
2151  const
2152 {
2153  //* Calculate distance to other particle [cm]
2154 
2155  Double_t dS, dS1;
2156  GetDStoParticle( p, dS, dS1 );
2157  Double_t mP[8], mC[36], mP1[8], mC1[36];
2158  Transport( dS, mP, mC );
2159  p.Transport( dS1, mP1, mC1 );
2160  Double_t dx = mP[0]-mP1[0];
2161  Double_t dy = mP[1]-mP1[1];
2162  Double_t dz = mP[2]-mP1[2];
2163  dz = 0;
2164  return TMath::Sqrt(dx*dx+dy*dy+dz*dz);
2165 }
2166 
2168 {
2169  //* Calculate sqrt(Chi2/ndf) deviation from vertex
2170 
2171  return GetDeviationFromVertex( Vtx.fP, Vtx.fC );
2172 }
2173 
2174 
2175 Double_t AliKFParticleBase::GetDeviationFromVertex( const Double_t v[], const Double_t Cv[] ) const
2176 {
2177  //* Calculate sqrt(Chi2/ndf) deviation from vertex
2178  //* v = [xyz], Cv=[Cxx,Cxy,Cyy,Cxz,Cyz,Czz]-covariance matrix
2179 
2180  Double_t mP[8];
2181  Double_t mC[36];
2182 
2183  Transport( GetDStoPoint(v), mP, mC );
2184 
2185  Double_t d[3]={ v[0]-mP[0], v[1]-mP[1], v[2]-mP[2]};
2186 
2187  Double_t sigmaS = .1+10.*TMath::Sqrt( (d[0]*d[0]+d[1]*d[1]+d[2]*d[2])/
2188  (mP[3]*mP[3]+mP[4]*mP[4]+mP[5]*mP[5]) );
2189 
2190 
2191  Double_t h[3] = { mP[3]*sigmaS, mP[4]*sigmaS, mP[5]*sigmaS };
2192 
2193  Double_t mSi[6] =
2194  { mC[0] +h[0]*h[0],
2195  mC[1] +h[1]*h[0], mC[2] +h[1]*h[1],
2196  mC[3] +h[2]*h[0], mC[4] +h[2]*h[1], mC[5] +h[2]*h[2] };
2197 
2198  if( Cv ){
2199  mSi[0]+=Cv[0];
2200  mSi[1]+=Cv[1];
2201  mSi[2]+=Cv[2];
2202  mSi[3]+=Cv[3];
2203  mSi[4]+=Cv[4];
2204  mSi[5]+=Cv[5];
2205  }
2206 
2207  Double_t mS[6];
2208 
2209  mS[0] = mSi[2]*mSi[5] - mSi[4]*mSi[4];
2210  mS[1] = mSi[3]*mSi[4] - mSi[1]*mSi[5];
2211  mS[2] = mSi[0]*mSi[5] - mSi[3]*mSi[3];
2212  mS[3] = mSi[1]*mSi[4] - mSi[2]*mSi[3];
2213  mS[4] = mSi[1]*mSi[3] - mSi[0]*mSi[4];
2214  mS[5] = mSi[0]*mSi[2] - mSi[1]*mSi[1];
2215 
2216  Double_t s = ( mSi[0]*mS[0] + mSi[1]*mS[1] + mSi[3]*mS[3] );
2217  s = ( s > 1.E-20 ) ?1./s :0;
2218 
2219  return TMath::Sqrt( TMath::Abs(s*( ( mS[0]*d[0] + mS[1]*d[1] + mS[3]*d[2])*d[0]
2220  +(mS[1]*d[0] + mS[2]*d[1] + mS[4]*d[2])*d[1]
2221  +(mS[3]*d[0] + mS[4]*d[1] + mS[5]*d[2])*d[2] ))/2);
2222 }
2223 
2224 
2226  const
2227 {
2228  //* Calculate sqrt(Chi2/ndf) deviation from other particle
2229 
2230  Double_t dS, dS1;
2231  GetDStoParticle( p, dS, dS1 );
2232  Double_t mP1[8], mC1[36];
2233  p.Transport( dS1, mP1, mC1 );
2234 
2235  Double_t d[3]={ fP[0]-mP1[0], fP[1]-mP1[1], fP[2]-mP1[2]};
2236 
2237  Double_t sigmaS = .1+10.*TMath::Sqrt( (d[0]*d[0]+d[1]*d[1]+d[2]*d[2])/
2238  (mP1[3]*mP1[3]+mP1[4]*mP1[4]+mP1[5]*mP1[5]) );
2239 
2240  Double_t h[3] = { mP1[3]*sigmaS, mP1[4]*sigmaS, mP1[5]*sigmaS };
2241 
2242  mC1[0] +=h[0]*h[0];
2243  mC1[1] +=h[1]*h[0];
2244  mC1[2] +=h[1]*h[1];
2245  mC1[3] +=h[2]*h[0];
2246  mC1[4] +=h[2]*h[1];
2247  mC1[5] +=h[2]*h[2];
2248 
2249  return GetDeviationFromVertex( mP1, mC1 )*TMath::Sqrt(2./1.);
2250 }
2251 
2252 
2253 
2255 {
2256  //* Subtract the particle from the vertex
2257 
2258  Double_t fld[3];
2259  {
2260  GetFieldValue( Vtx.fP, fld );
2261  const Double_t kCLight = 0.000299792458;
2262  fld[0]*=kCLight; fld[1]*=kCLight; fld[2]*=kCLight;
2263  }
2264 
2265  Double_t m[8];
2266  Double_t mCm[36];
2267 
2268  if( Vtx.fIsLinearized ){
2269  GetMeasurement( Vtx.fVtxGuess, m, mCm );
2270  } else {
2271  GetMeasurement( Vtx.fP, m, mCm );
2272  }
2273 
2274  Double_t mV[6];
2275 
2276  mV[ 0] = mCm[ 0];
2277  mV[ 1] = mCm[ 1];
2278  mV[ 2] = mCm[ 2];
2279  mV[ 3] = mCm[ 3];
2280  mV[ 4] = mCm[ 4];
2281  mV[ 5] = mCm[ 5];
2282 
2283  //*
2284 
2285  Double_t mS[6];
2286  {
2287  Double_t mSi[6] = { mV[0]-Vtx.fC[0],
2288  mV[1]-Vtx.fC[1], mV[2]-Vtx.fC[2],
2289  mV[3]-Vtx.fC[3], mV[4]-Vtx.fC[4], mV[5]-Vtx.fC[5] };
2290 
2291  mS[0] = mSi[2]*mSi[5] - mSi[4]*mSi[4];
2292  mS[1] = mSi[3]*mSi[4] - mSi[1]*mSi[5];
2293  mS[2] = mSi[0]*mSi[5] - mSi[3]*mSi[3];
2294  mS[3] = mSi[1]*mSi[4] - mSi[2]*mSi[3];
2295  mS[4] = mSi[1]*mSi[3] - mSi[0]*mSi[4];
2296  mS[5] = mSi[0]*mSi[2] - mSi[1]*mSi[1];
2297 
2298  Double_t s = ( mSi[0]*mS[0] + mSi[1]*mS[1] + mSi[3]*mS[3] );
2299  s = ( s > 1.E-20 ) ?1./s :0;
2300  mS[0]*=s;
2301  mS[1]*=s;
2302  mS[2]*=s;
2303  mS[3]*=s;
2304  mS[4]*=s;
2305  mS[5]*=s;
2306  }
2307 
2308  //* Residual (measured - estimated)
2309 
2310  Double_t zeta[3] = { m[0]-Vtx.fP[0], m[1]-Vtx.fP[1], m[2]-Vtx.fP[2] };
2311 
2312  //* mCHt = mCH' - D'
2313 
2314  Double_t mCHt0[3], mCHt1[3], mCHt2[3];
2315 
2316  mCHt0[0]=Vtx.fC[ 0] ; mCHt1[0]=Vtx.fC[ 1] ; mCHt2[0]=Vtx.fC[ 3] ;
2317  mCHt0[1]=Vtx.fC[ 1] ; mCHt1[1]=Vtx.fC[ 2] ; mCHt2[1]=Vtx.fC[ 4] ;
2318  mCHt0[2]=Vtx.fC[ 3] ; mCHt1[2]=Vtx.fC[ 4] ; mCHt2[2]=Vtx.fC[ 5] ;
2319 
2320  //* Kalman gain K = mCH'*S
2321 
2322  Double_t k0[3], k1[3], k2[3];
2323 
2324  for(Int_t i=0;i<3;++i){
2325  k0[i] = mCHt0[i]*mS[0] + mCHt1[i]*mS[1] + mCHt2[i]*mS[3];
2326  k1[i] = mCHt0[i]*mS[1] + mCHt1[i]*mS[2] + mCHt2[i]*mS[4];
2327  k2[i] = mCHt0[i]*mS[3] + mCHt1[i]*mS[4] + mCHt2[i]*mS[5];
2328  }
2329 
2330  //* New estimation of the vertex position r += K*zeta
2331 
2332  Double_t dChi2 = -(mS[0]*zeta[0] + mS[1]*zeta[1] + mS[3]*zeta[2])*zeta[0]
2333  + (mS[1]*zeta[0] + mS[2]*zeta[1] + mS[4]*zeta[2])*zeta[1]
2334  + (mS[3]*zeta[0] + mS[4]*zeta[1] + mS[5]*zeta[2])*zeta[2];
2335 
2336  if( Vtx.fChi2 - dChi2 < 0 ) return;
2337 
2338  for(Int_t i=0;i<3;++i)
2339  Vtx.fP[i] -= k0[i]*zeta[0] + k1[i]*zeta[1] + k2[i]*zeta[2];
2340 
2341  //* New covariance matrix C -= K*(mCH')'
2342 
2343  for(Int_t i=0, k=0;i<3;++i){
2344  for(Int_t j=0;j<=i;++j,++k)
2345  Vtx.fC[k] += k0[i]*mCHt0[j] + k1[i]*mCHt1[j] + k2[i]*mCHt2[j];
2346  }
2347 
2348  //* Calculate Chi^2
2349 
2350  Vtx.fNDF -= 2;
2351  Vtx.fChi2 -= dChi2;
2352 }
2353 
2354 
2355 
2357  Double_t P[], Double_t C[] ) const
2358 {
2359  //* Transport the particle as a straight line
2360 
2361  P[0] = fP[0] + dS*fP[3];
2362  P[1] = fP[1] + dS*fP[4];
2363  P[2] = fP[2] + dS*fP[5];
2364  P[3] = fP[3];
2365  P[4] = fP[4];
2366  P[5] = fP[5];
2367  P[6] = fP[6];
2368  P[7] = fP[7];
2369 
2370  Double_t c6 = fC[ 6] + dS*fC[ 9];
2371  Double_t c11 = fC[11] + dS*fC[14];
2372  Double_t c17 = fC[17] + dS*fC[20];
2373  Double_t sc13 = dS*fC[13];
2374  Double_t sc18 = dS*fC[18];
2375  Double_t sc19 = dS*fC[19];
2376 
2377  C[ 0] = fC[ 0] + dS*( fC[ 6] + c6 );
2378  C[ 2] = fC[ 2] + dS*( fC[11] + c11 );
2379  C[ 5] = fC[ 5] + dS*( fC[17] + c17 );
2380 
2381  C[ 7] = fC[ 7] + sc13;
2382  C[ 8] = fC[ 8] + sc18;
2383  C[ 9] = fC[ 9];
2384 
2385  C[12] = fC[12] + sc19;
2386 
2387  C[ 1] = fC[ 1] + dS*( fC[10] + C[ 7] );
2388  C[ 3] = fC[ 3] + dS*( fC[15] + C[ 8] );
2389  C[ 4] = fC[ 4] + dS*( fC[16] + C[12] );
2390  C[ 6] = c6;
2391 
2392  C[10] = fC[10] + sc13;
2393  C[11] = c11;
2394 
2395  C[13] = fC[13];
2396  C[14] = fC[14];
2397  C[15] = fC[15] + sc18;
2398  C[16] = fC[16] + sc19;
2399  C[17] = c17;
2400 
2401  C[18] = fC[18];
2402  C[19] = fC[19];
2403  C[20] = fC[20];
2404  C[21] = fC[21] + dS*fC[24];
2405  C[22] = fC[22] + dS*fC[25];
2406  C[23] = fC[23] + dS*fC[26];
2407 
2408  C[24] = fC[24];
2409  C[25] = fC[25];
2410  C[26] = fC[26];
2411  C[27] = fC[27];
2412  C[28] = fC[28] + dS*fC[31];
2413  C[29] = fC[29] + dS*fC[32];
2414  C[30] = fC[30] + dS*fC[33];
2415 
2416  C[31] = fC[31];
2417  C[32] = fC[32];
2418  C[33] = fC[33];
2419  C[34] = fC[34];
2420  C[35] = fC[35];
2421 }
2422 
2423 
2425  const AliKFParticleBase &daughter2, double Bz )
2426 {
2427  //* Create gamma
2428 
2429  const AliKFParticleBase *daughters[2] = { &daughter1, &daughter2};
2430 
2431  double v0[3];
2432 
2433  if( !fIsLinearized ){
2434  Double_t ds, ds1;
2435  Double_t m[8];
2436  Double_t mCd[36];
2437  daughter1.GetDStoParticle(daughter2, ds, ds1);
2438  daughter1.Transport( ds, m, mCd );
2439  fP[0] = m[0];
2440  fP[1] = m[1];
2441  fP[2] = m[2];
2442  daughter2.Transport( ds1, m, mCd );
2443  fP[0] = .5*( fP[0] + m[0] );
2444  fP[1] = .5*( fP[1] + m[1] );
2445  fP[2] = .5*( fP[2] + m[2] );
2446  } else {
2447  fP[0] = fVtxGuess[0];
2448  fP[1] = fVtxGuess[1];
2449  fP[2] = fVtxGuess[2];
2450  }
2451 
2452  double daughterP[2][8], daughterC[2][36];
2453  double vtxMom[2][3];
2454 
2455  int nIter = fIsLinearized ?1 :2;
2456 
2457  for( int iter=0; iter<nIter; iter++){
2458 
2459  v0[0] = fP[0];
2460  v0[1] = fP[1];
2461  v0[2] = fP[2];
2462 
2463  fAtProductionVertex = 0;
2464  fSFromDecay = 0;
2465  fP[0] = v0[0];
2466  fP[1] = v0[1];
2467  fP[2] = v0[2];
2468  fP[3] = 0;
2469  fP[4] = 0;
2470  fP[5] = 0;
2471  fP[6] = 0;
2472  fP[7] = 0;
2473 
2474 
2475  // fit daughters to the vertex guess
2476 
2477  {
2478  for( int id=0; id<2; id++ ){
2479 
2480  double *p = daughterP[id];
2481  double *mC = daughterC[id];
2482 
2483  daughters[id]->GetMeasurement( v0, p, mC );
2484 
2485  Double_t mAi[6];
2486  InvertSym3(mC, mAi );
2487 
2488  Double_t mB[3][3];
2489 
2490  mB[0][0] = mC[ 6]*mAi[0] + mC[ 7]*mAi[1] + mC[ 8]*mAi[3];
2491  mB[0][1] = mC[ 6]*mAi[1] + mC[ 7]*mAi[2] + mC[ 8]*mAi[4];
2492  mB[0][2] = mC[ 6]*mAi[3] + mC[ 7]*mAi[4] + mC[ 8]*mAi[5];
2493 
2494  mB[1][0] = mC[10]*mAi[0] + mC[11]*mAi[1] + mC[12]*mAi[3];
2495  mB[1][1] = mC[10]*mAi[1] + mC[11]*mAi[2] + mC[12]*mAi[4];
2496  mB[1][2] = mC[10]*mAi[3] + mC[11]*mAi[4] + mC[12]*mAi[5];
2497 
2498  mB[2][0] = mC[15]*mAi[0] + mC[16]*mAi[1] + mC[17]*mAi[3];
2499  mB[2][1] = mC[15]*mAi[1] + mC[16]*mAi[2] + mC[17]*mAi[4];
2500  mB[2][2] = mC[15]*mAi[3] + mC[16]*mAi[4] + mC[17]*mAi[5];
2501 
2502  Double_t z[3] = { v0[0]-p[0], v0[1]-p[1], v0[2]-p[2] };
2503 
2504  vtxMom[id][0] = p[3] + mB[0][0]*z[0] + mB[0][1]*z[1] + mB[0][2]*z[2];
2505  vtxMom[id][1] = p[4] + mB[1][0]*z[0] + mB[1][1]*z[1] + mB[1][2]*z[2];
2506  vtxMom[id][2] = p[5] + mB[2][0]*z[0] + mB[2][1]*z[1] + mB[2][2]*z[2];
2507 
2508  daughters[id]->Transport( daughters[id]->GetDStoPoint(v0), p, mC );
2509 
2510  }
2511 
2512  } // fit daughters to guess
2513 
2514 
2515  // fit new vertex
2516  {
2517 
2518  double mpx0 = vtxMom[0][0]+vtxMom[1][0];
2519  double mpy0 = vtxMom[0][1]+vtxMom[1][1];
2520  double mpt0 = TMath::Sqrt(mpx0*mpx0 + mpy0*mpy0);
2521  // double a0 = TMath::ATan2(mpy0,mpx0);
2522 
2523  double ca0 = mpx0/mpt0;
2524  double sa0 = mpy0/mpt0;
2525  double r[3] = { v0[0], v0[1], v0[2] };
2526  double mC[3][3] = {{1000., 0 , 0 },
2527  {0, 1000., 0 },
2528  {0, 0, 1000. } };
2529  double chi2=0;
2530 
2531  for( int id=0; id<2; id++ ){
2532  const Double_t kCLight = 0.000299792458;
2533  Double_t q = Bz*daughters[id]->GetQ()*kCLight;
2534  Double_t px0 = vtxMom[id][0];
2535  Double_t py0 = vtxMom[id][1];
2536  Double_t pz0 = vtxMom[id][2];
2537  Double_t pt0 = TMath::Sqrt(px0*px0+py0*py0);
2538  Double_t mG[3][6], mB[3], mH[3][3];
2539  // r = {vx,vy,vz};
2540  // m = {x,y,z,Px,Py,Pz};
2541  // V = daughter.C
2542  // G*m + B = H*r;
2543  // q*x + Py - q*vx - sin(a)*Pt = 0
2544  // q*y - Px - q*vy + cos(a)*Pt = 0
2545  // (Px*cos(a) + Py*sin(a) ) (vz -z) - Pz( cos(a)*(vx-x) + sin(a)*(vy-y)) = 0
2546 
2547  mG[0][0] = q;
2548  mG[0][1] = 0;
2549  mG[0][2] = 0;
2550  mG[0][3] = -sa0*px0/pt0;
2551  mG[0][4] = 1 -sa0*py0/pt0;
2552  mG[0][5] = 0;
2553  mH[0][0] = q;
2554  mH[0][1] = 0;
2555  mH[0][2] = 0;
2556  mB[0] = py0 - sa0*pt0 - mG[0][3]*px0 - mG[0][4]*py0 ;
2557 
2558  // q*y - Px - q*vy + cos(a)*Pt = 0
2559 
2560  mG[1][0] = 0;
2561  mG[1][1] = q;
2562  mG[1][2] = 0;
2563  mG[1][3] = -1 + ca0*px0/pt0;
2564  mG[1][4] = + ca0*py0/pt0;
2565  mG[1][5] = 0;
2566  mH[1][0] = 0;
2567  mH[1][1] = q;
2568  mH[1][2] = 0;
2569  mB[1] = -px0 + ca0*pt0 - mG[1][3]*px0 - mG[1][4]*py0 ;
2570 
2571  // (Px*cos(a) + Py*sin(a) ) (z -vz) - Pz( cos(a)*(x-vx) + sin(a)*(y-vy)) = 0
2572 
2573  mG[2][0] = -pz0*ca0;
2574  mG[2][1] = -pz0*sa0;
2575  mG[2][2] = px0*ca0 + py0*sa0;
2576  mG[2][3] = 0;
2577  mG[2][4] = 0;
2578  mG[2][5] = 0;
2579 
2580  mH[2][0] = mG[2][0];
2581  mH[2][1] = mG[2][1];
2582  mH[2][2] = mG[2][2];
2583 
2584  mB[2] = 0;
2585 
2586  // fit the vertex
2587 
2588  // V = GVGt
2589 
2590  double mGV[3][6];
2591  double mV[6];
2592  double m[3];
2593  for( int i=0; i<3; i++ ){
2594  m[i] = mB[i];
2595  for( int k=0; k<6; k++ ) m[i]+=mG[i][k]*daughterP[id][k];
2596  }
2597  for( int i=0; i<3; i++ ){
2598  for( int j=0; j<6; j++ ){
2599  mGV[i][j] = 0;
2600  for( int k=0; k<6; k++ ) mGV[i][j]+=mG[i][k]*daughterC[id][ IJ(k,j) ];
2601  }
2602  }
2603  for( int i=0, k=0; i<3; i++ ){
2604  for( int j=0; j<=i; j++,k++ ){
2605  mV[k] = 0;
2606  for( int l=0; l<6; l++ ) mV[k]+=mGV[i][l]*mG[j][l];
2607  }
2608  }
2609 
2610 
2611  //* CHt
2612 
2613  Double_t mCHt[3][3];
2614  Double_t mHCHt[6];
2615  Double_t mHr[3];
2616  for( int i=0; i<3; i++ ){
2617  mHr[i] = 0;
2618  for( int k=0; k<3; k++ ) mHr[i]+= mH[i][k]*r[k];
2619  }
2620 
2621  for( int i=0; i<3; i++ ){
2622  for( int j=0; j<3; j++){
2623  mCHt[i][j] = 0;
2624  for( int k=0; k<3; k++ ) mCHt[i][j]+= mC[i][k]*mH[j][k];
2625  }
2626  }
2627 
2628  for( int i=0, k=0; i<3; i++ ){
2629  for( int j=0; j<=i; j++, k++ ){
2630  mHCHt[k] = 0;
2631  for( int l=0; l<3; l++ ) mHCHt[k]+= mH[i][l]*mCHt[l][j];
2632  }
2633  }
2634 
2635  Double_t mS[6] = { mHCHt[0]+mV[0],
2636  mHCHt[1]+mV[1], mHCHt[2]+mV[2],
2637  mHCHt[3]+mV[3], mHCHt[4]+mV[4], mHCHt[5]+mV[5] };
2638 
2639 
2640  InvertSym3(mS,mS);
2641 
2642  //* Residual (measured - estimated)
2643 
2644  Double_t zeta[3] = { m[0]-mHr[0], m[1]-mHr[1], m[2]-mHr[2] };
2645 
2646  //* Kalman gain K = mCH'*S
2647 
2648  Double_t k[3][3];
2649 
2650  for(Int_t i=0;i<3;++i){
2651  k[i][0] = mCHt[i][0]*mS[0] + mCHt[i][1]*mS[1] + mCHt[i][2]*mS[3];
2652  k[i][1] = mCHt[i][0]*mS[1] + mCHt[i][1]*mS[2] + mCHt[i][2]*mS[4];
2653  k[i][2] = mCHt[i][0]*mS[3] + mCHt[i][1]*mS[4] + mCHt[i][2]*mS[5];
2654  }
2655 
2656  //* New estimation of the vertex position r += K*zeta
2657 
2658  for(Int_t i=0;i<3;++i)
2659  r[i] = r[i] + k[i][0]*zeta[0] + k[i][1]*zeta[1] + k[i][2]*zeta[2];
2660 
2661  //* New covariance matrix C -= K*(mCH')'
2662 
2663  for(Int_t i=0;i<3;++i){
2664  for(Int_t j=0;j<=i;++j){
2665  mC[i][j] = mC[i][j] - (k[i][0]*mCHt[j][0] + k[i][1]*mCHt[j][1] + k[i][2]*mCHt[j][2]);
2666  mC[j][i] = mC[i][j];
2667  }
2668  }
2669 
2670  //* Calculate Chi^2
2671 
2672  chi2 += ( ( mS[0]*zeta[0] + mS[1]*zeta[1] + mS[3]*zeta[2] )*zeta[0]
2673  +(mS[1]*zeta[0] + mS[2]*zeta[1] + mS[4]*zeta[2] )*zeta[1]
2674  +(mS[3]*zeta[0] + mS[4]*zeta[1] + mS[5]*zeta[2] )*zeta[2] );
2675  }
2676 
2677  // store vertex
2678 
2679  fNDF = 2;
2680  fChi2 = chi2;
2681  for( int i=0; i<3; i++ ) fP[i] = r[i];
2682  for( int i=0,k=0; i<3; i++ ){
2683  for( int j=0; j<=i; j++,k++ ){
2684  fC[k] = mC[i][j];
2685  }
2686  }
2687  }
2688 
2689  } // iterations
2690 
2691  // now fit daughters to the vertex
2692 
2693  fQ = 0;
2694  fSFromDecay = 0;
2695 
2696  for(Int_t i=3;i<8;++i) fP[i]=0.;
2697  for(Int_t i=6;i<35;++i) fC[i]=0.;
2698  fC[35] = 100.;
2699 
2700  for( int id=0; id<2; id++ ){
2701 
2702  double *p = daughterP[id];
2703  double *mC = daughterC[id];
2704  daughters[id]->GetMeasurement( v0, p, mC );
2705 
2706  const Double_t *m = fP, *mV = fC;
2707 
2708  Double_t mAi[6];
2709  InvertSym3(mC, mAi );
2710 
2711  Double_t mB[4][3];
2712 
2713  mB[0][0] = mC[ 6]*mAi[0] + mC[ 7]*mAi[1] + mC[ 8]*mAi[3];
2714  mB[0][1] = mC[ 6]*mAi[1] + mC[ 7]*mAi[2] + mC[ 8]*mAi[4];
2715  mB[0][2] = mC[ 6]*mAi[3] + mC[ 7]*mAi[4] + mC[ 8]*mAi[5];
2716 
2717  mB[1][0] = mC[10]*mAi[0] + mC[11]*mAi[1] + mC[12]*mAi[3];
2718  mB[1][1] = mC[10]*mAi[1] + mC[11]*mAi[2] + mC[12]*mAi[4];
2719  mB[1][2] = mC[10]*mAi[3] + mC[11]*mAi[4] + mC[12]*mAi[5];
2720 
2721  mB[2][0] = mC[15]*mAi[0] + mC[16]*mAi[1] + mC[17]*mAi[3];
2722  mB[2][1] = mC[15]*mAi[1] + mC[16]*mAi[2] + mC[17]*mAi[4];
2723  mB[2][2] = mC[15]*mAi[3] + mC[16]*mAi[4] + mC[17]*mAi[5];
2724 
2725  mB[3][0] = mC[21]*mAi[0] + mC[22]*mAi[1] + mC[23]*mAi[3];
2726  mB[3][1] = mC[21]*mAi[1] + mC[22]*mAi[2] + mC[23]*mAi[4];
2727  mB[3][2] = mC[21]*mAi[3] + mC[22]*mAi[4] + mC[23]*mAi[5];
2728 
2729 
2730  Double_t z[3] = { m[0]-p[0], m[1]-p[1], m[2]-p[2] };
2731 
2732 // {
2733 // Double_t mAV[6] = { mC[0]-mV[0], mC[1]-mV[1], mC[2]-mV[2],
2734 // mC[3]-mV[3], mC[4]-mV[4], mC[5]-mV[5] };
2735 //
2736 // Double_t mAVi[6];
2737 // if( !InvertSym3(mAV, mAVi) ){
2738 // Double_t dChi2 = ( +(mAVi[0]*z[0] + mAVi[1]*z[1] + mAVi[3]*z[2])*z[0]
2739 // +(mAVi[1]*z[0] + mAVi[2]*z[1] + mAVi[4]*z[2])*z[1]
2740 // +(mAVi[3]*z[0] + mAVi[4]*z[1] + mAVi[5]*z[2])*z[2] );
2741 // fChi2+= TMath::Abs( dChi2 );
2742 // }
2743 // fNDF += 2;
2744 // }
2745 
2746  //* Add the daughter momentum to the particle momentum
2747 
2748  fP[3]+= p[3] + mB[0][0]*z[0] + mB[0][1]*z[1] + mB[0][2]*z[2];
2749  fP[4]+= p[4] + mB[1][0]*z[0] + mB[1][1]*z[1] + mB[1][2]*z[2];
2750  fP[5]+= p[5] + mB[2][0]*z[0] + mB[2][1]*z[1] + mB[2][2]*z[2];
2751  fP[6]+= p[6] + mB[3][0]*z[0] + mB[3][1]*z[1] + mB[3][2]*z[2];
2752 
2753  Double_t d0, d1, d2;
2754 
2755  d0= mB[0][0]*mV[0] + mB[0][1]*mV[1] + mB[0][2]*mV[3] - mC[ 6];
2756  d1= mB[0][0]*mV[1] + mB[0][1]*mV[2] + mB[0][2]*mV[4] - mC[ 7];
2757  d2= mB[0][0]*mV[3] + mB[0][1]*mV[4] + mB[0][2]*mV[5] - mC[ 8];
2758 
2759  //fC[6]+= mC[ 6] + d0;
2760  //fC[7]+= mC[ 7] + d1;
2761  //fC[8]+= mC[ 8] + d2;
2762  fC[9]+= mC[ 9] + d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
2763 
2764  d0= mB[1][0]*mV[0] + mB[1][1]*mV[1] + mB[1][2]*mV[3] - mC[10];
2765  d1= mB[1][0]*mV[1] + mB[1][1]*mV[2] + mB[1][2]*mV[4] - mC[11];
2766  d2= mB[1][0]*mV[3] + mB[1][1]*mV[4] + mB[1][2]*mV[5] - mC[12];
2767 
2768  //fC[10]+= mC[10]+ d0;
2769  //fC[11]+= mC[11]+ d1;
2770  //fC[12]+= mC[12]+ d2;
2771  fC[13]+= mC[13]+ d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
2772  fC[14]+= mC[14]+ d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
2773 
2774  d0= mB[2][0]*mV[0] + mB[2][1]*mV[1] + mB[2][2]*mV[3] - mC[15];
2775  d1= mB[2][0]*mV[1] + mB[2][1]*mV[2] + mB[2][2]*mV[4] - mC[16];
2776  d2= mB[2][0]*mV[3] + mB[2][1]*mV[4] + mB[2][2]*mV[5] - mC[17];
2777 
2778  //fC[15]+= mC[15]+ d0;
2779  //fC[16]+= mC[16]+ d1;
2780  //fC[17]+= mC[17]+ d2;
2781  fC[18]+= mC[18]+ d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
2782  fC[19]+= mC[19]+ d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
2783  fC[20]+= mC[20]+ d0*mB[2][0] + d1*mB[2][1] + d2*mB[2][2];
2784 
2785  d0= mB[3][0]*mV[0] + mB[3][1]*mV[1] + mB[3][2]*mV[3] - mC[21];
2786  d1= mB[3][0]*mV[1] + mB[3][1]*mV[2] + mB[3][2]*mV[4] - mC[22];
2787  d2= mB[3][0]*mV[3] + mB[3][1]*mV[4] + mB[3][2]*mV[5] - mC[23];
2788 
2789  //fC[21]+= mC[21] + d0;
2790  //fC[22]+= mC[22] + d1;
2791  //fC[23]+= mC[23] + d2;
2792  fC[24]+= mC[24] + d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
2793  fC[25]+= mC[25] + d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
2794  fC[26]+= mC[26] + d0*mB[2][0] + d1*mB[2][1] + d2*mB[2][2];
2795  fC[27]+= mC[27] + d0*mB[3][0] + d1*mB[3][1] + d2*mB[3][2];
2796  }
2797 
2798 // SetMassConstraint(0,0);
2800 }
2801 
2803 {
2804 // example:
2805 // AliKFParticle PosParticle(...)
2806 // AliKFParticle NegParticle(...)
2807 // Gamma.ConstructGamma(PosParticle, NegParticle);
2808 // Double_t VertexGamma[3] = {Gamma.GetX(), Gamma.GetY(), Gamma.GetZ()};
2809 // PosParticle.TransportToPoint(VertexGamma);
2810 // NegParticle.TransportToPoint(VertexGamma);
2811 // Double_t armenterosQtAlfa[2] = {0.};
2812 // AliKFParticle::GetArmenterosPodolanski(PosParticle, NegParticle, armenterosQtAlfa );
2813 
2814  Double_t alpha = 0., qt = 0.;
2815  Double_t spx = positive.GetPx() + negative.GetPx();
2816  Double_t spy = positive.GetPy() + negative.GetPy();
2817  Double_t spz = positive.GetPz() + negative.GetPz();
2818  Double_t sp = sqrt(spx*spx + spy*spy + spz*spz);
2819  if( sp == 0.0) return;
2820  Double_t pn, pln, plp; // ,pp;
2821 
2822  pn = TMath::Sqrt(negative.GetPx()*negative.GetPx() + negative.GetPy()*negative.GetPy() + negative.GetPz()*negative.GetPz());
2823  // pp = TMath::Sqrt(positive.GetPx()*positive.GetPx() + positive.GetPy()*positive.GetPy() + positive.GetPz()*positive.GetPz());
2824  pln = (negative.GetPx()*spx+negative.GetPy()*spy+negative.GetPz()*spz)/sp;
2825  plp = (positive.GetPx()*spx+positive.GetPy()*spy+positive.GetPz()*spz)/sp;
2826 
2827  if( pn == 0.0) return;
2828  Double_t ptm = (1.-((pln/pn)*(pln/pn)));
2829  qt= (ptm>=0.)? pn*sqrt(ptm) :0;
2830  alpha = (plp-pln)/(plp+pln);
2831 
2832  QtAlfa[0] = qt;
2833  QtAlfa[1] = alpha;
2834 }
2835 
2836 void AliKFParticleBase::RotateXY(Double_t angle, Double_t Vtx[3])
2837 {
2838  // Rotates the KFParticle object around OZ axis, OZ axis is set by the vertex position
2839  // Double_t angle - angle of rotation in XY plane in [rad]
2840  // Double_t Vtx[3] - position of the vertex in [cm]
2841 
2842  // Before rotation the center of the coordinat system should be moved to the vertex position; move back after rotation
2843  X() = X() - Vtx[0];
2844  Y() = Y() - Vtx[1];
2845  Z() = Z() - Vtx[2];
2846 
2847  // Rotate the kf particle
2848  Double_t c = TMath::Cos(angle);
2849  Double_t s = TMath::Sin(angle);
2850 
2851  Double_t mA[8][ 8];
2852  for( Int_t i=0; i<8; i++ ){
2853  for( Int_t j=0; j<8; j++){
2854  mA[i][j] = 0;
2855  }
2856  }
2857  for( int i=0; i<8; i++ ){
2858  mA[i][i] = 1;
2859  }
2860  mA[0][0] = c; mA[0][1] = s;
2861  mA[1][0] = -s; mA[1][1] = c;
2862  mA[3][3] = c; mA[3][4] = s;
2863  mA[4][3] = -s; mA[4][4] = c;
2864 
2865  Double_t mAC[8][8];
2866  Double_t mAp[8];
2867 
2868  for( Int_t i=0; i<8; i++ ){
2869  mAp[i] = 0;
2870  for( Int_t k=0; k<8; k++){
2871  mAp[i]+=mA[i][k] * fP[k];
2872  }
2873  }
2874 
2875  for( Int_t i=0; i<8; i++){
2876  fP[i] = mAp[i];
2877  }
2878 
2879  for( Int_t i=0; i<8; i++ ){
2880  for( Int_t j=0; j<8; j++ ){
2881  mAC[i][j] = 0;
2882  for( Int_t k=0; k<8; k++ ){
2883  mAC[i][j]+= mA[i][k] * GetCovariance(k,j);
2884  }
2885  }
2886  }
2887 
2888  for( Int_t i=0; i<8; i++ ){
2889  for( Int_t j=0; j<=i; j++ ){
2890  Double_t xx = 0;
2891  for( Int_t k=0; k<8; k++){
2892  xx+= mAC[i][k]*mA[j][k];
2893  }
2894  Covariance(i,j) = xx;
2895  }
2896  }
2897 
2898  X() = GetX() + Vtx[0];
2899  Y() = GetY() + Vtx[1];
2900  Z() = GetZ() + Vtx[2];
2901 }
2902 
2903 Bool_t AliKFParticleBase::InvertSym3( const Double_t A[], Double_t Ai[] )
2904 {
2905  //* Invert symmetric matric stored in low-triagonal form
2906 
2907  bool ret = 0;
2908  double a0 = A[0], a1 = A[1], a2 = A[2], a3 = A[3];
2909 
2910  Ai[0] = a2*A[5] - A[4]*A[4];
2911  Ai[1] = a3*A[4] - a1*A[5];
2912  Ai[3] = a1*A[4] - a2*a3;
2913  Double_t det = (a0*Ai[0] + a1*Ai[1] + a3*Ai[3]);
2914  if( TMath::Abs(det)>1.e-20 ) det = 1./det;
2915  else{
2916  det = 0;
2917  ret = 1;
2918  }
2919  Ai[0] *= det;
2920  Ai[1] *= det;
2921  Ai[3] *= det;
2922  Ai[2] = ( a0*A[5] - a3*a3 )*det;
2923  Ai[4] = ( a1*a3 - a0*A[4] )*det;
2924  Ai[5] = ( a0*a2 - a1*a1 )*det;
2925  return ret;
2926 }
2927 
2928 void AliKFParticleBase::MultQSQt( const Double_t Q[], const Double_t S[], Double_t SOut[] )
2929 {
2930  //* Matrix multiplication Q*S*Q^T, Q - square matrix, S - symmetric
2931 
2932  const Int_t kN= 8;
2933  Double_t mA[kN*kN];
2934 
2935  for( Int_t i=0, ij=0; i<kN; i++ ){
2936  for( Int_t j=0; j<kN; j++, ++ij ){
2937  mA[ij] = 0 ;
2938  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];
2939  }
2940  }
2941 
2942  for( Int_t i=0; i<kN; i++ ){
2943  for( Int_t j=0; j<=i; j++ ){
2944  Int_t ij = ( j<=i ) ? i*(i+1)/2+j :j*(j+1)/2+i;
2945  SOut[ij] = 0 ;
2946  for( Int_t k=0; k<kN; k++ ) SOut[ij] += Q[ i*kN+k ] * mA[ k*kN+j ];
2947  }
2948  }
2949 }
2950 
2951 
2952 // 72-charachters line to define the printer border
2953 //3456789012345678901234567890123456789012345678901234567890123456789012
2954 
TBrowser b
Definition: RunAnaESD.C:12
Int_t GetQ() const
Int_t GetPt(Double_t &Pt, Double_t &SigmaPt) const
static Bool_t InvertSym3(const Double_t A[], Double_t Ainv[])
void TransportLine(Double_t S, Double_t P[], Double_t C[]) const
static Int_t IJ(Int_t i, Int_t j)
Double_t GetDeviationFromParticle(const AliKFParticleBase &p) const
void AddDaughterWithEnergyFitMC(const AliKFParticleBase &Daughter)
Int_t GetPhi(Double_t &Phi, Double_t &SigmaPhi) const
Double_t GetDStoPointLine(const Double_t xyz[]) const
void RotateXY(Double_t angle, Double_t Vtx[3])
Double_t GetPz() const
const char * df
Definition: Raw2ESD.C:14
const Double_t & Y() const
void Construct(const AliKFParticleBase *vDaughters[], Int_t NDaughters, const AliKFParticleBase *ProdVtx=0, Double_t Mass=-1, Bool_t IsConstrained=0)
static Double_t GetSCorrection(const Double_t Part[], const Double_t XYZ[])
Double_t & Cij(Int_t i, Int_t j)
Float_t p[]
Definition: kNNTest.C:133
virtual Double_t GetDStoPoint(const Double_t xyz[]) const =0
Double_t GetY() const
void operator+=(const AliKFParticleBase &Daughter)
const Double_t & X() const
void AddDaughterWithEnergyCalc(const AliKFParticleBase &Daughter)
void ConstructGammaBz(const AliKFParticleBase &daughter1, const AliKFParticleBase &daughter2, double Bz)
Double_t GetDStoPointBz(Double_t Bz, const Double_t xyz[]) const
void SetProductionVertex(const AliKFParticleBase &Vtx)
Int_t GetDecayLength(Double_t &L, Double_t &SigmaL) const
Int_t GetLifeTime(Double_t &T, Double_t &SigmaT) const
const Int_t & Q() const
const Double_t & Z() const
void SetMassConstraint(Double_t Mass, Double_t SigmaMass=0)
Double_t chi2
Definition: AnalyzeLaser.C:7
Double_t GetDistanceFromVertex(const Double_t vtx[]) const
Double_t GetCovariance(Int_t i) const
Double_t GetX() const
Int_t GetEta(Double_t &Eta, Double_t &SigmaEta) const
void SubtractFromVertex(AliKFParticleBase &Vtx) const
virtual void GetDStoParticle(const AliKFParticleBase &p, Double_t &DS, Double_t &DSp) const =0
Double_t GetPy() const
void SetNonlinearMassConstraint(Double_t Mass)
void GetMeasurement(const Double_t XYZ[], Double_t m[], Double_t V[]) const
static void GetArmenterosPodolanski(AliKFParticleBase &positive, AliKFParticleBase &negative, Double_t QtAlfa[2])
Double_t GetDeviationFromVertex(const Double_t v[], const Double_t Cv[]=0) const
void TransportCBM(Double_t dS, Double_t P[], Double_t C[]) const
void TransportBz(Double_t Bz, Double_t dS, Double_t P[], Double_t C[]) const
const Double_t & S() const
virtual void Transport(Double_t dS, Double_t P[], Double_t C[]) const =0
Double_t GetZ() const
Int_t GetR(Double_t &R, Double_t &SigmaR) const
Double_t & Covariance(Int_t i)
AliTPCcalibV0 v0
void AddDaughterWithEnergyFit(const AliKFParticleBase &Daughter)
TF1 * f
Definition: interpolTest.C:21
virtual void GetFieldValue(const Double_t xyz[], Double_t B[]) const =0
void Convert(bool ToProduction)
Int_t GetMomentum(Double_t &P, Double_t &SigmaP) const
void SetVtxGuess(Double_t x, Double_t y, Double_t z)
void AddDaughter(const AliKFParticleBase &Daughter)
Double_t GetPx() const
Double_t GetDistanceFromParticle(const AliKFParticleBase &p) const
Int_t GetDecayLengthXY(Double_t &L, Double_t &SigmaL) const
static void MultQSQt(const Double_t Q[], const Double_t S[], Double_t SOut[])
void GetDStoParticleBz(Double_t Bz, const AliKFParticleBase &p, Double_t &dS, Double_t &dS1) const
void TransportToDS(Double_t dS)
Int_t GetMass(Double_t &M, Double_t &SigmaM) const