BmnRoot
Loading...
Searching...
No Matches
KFParticleBaseSIMD.cxx
Go to the documentation of this file.
1//---------------------------------------------------------------------------------
2// Implementation of the KFParticleBaseSIMD 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 KFParticle class
14//
15// -= Copyright &copy ALICE HLT and CBM L1 Groups =-
16//_________________________________________________________________________________
17
18
19#include "KFParticleBaseSIMD.h"
20
21#include <iostream>
22
23static const fvec Zero = 0;
24static const fvec One = 1;
25static const fvec small = 1.e-20;
26
27KFParticleBaseSIMD::KFParticleBaseSIMD() :fQ(0), fNDF(-3), fChi2(0), fSFromDecay(0), fAtProductionVertex(0), fIsVtxGuess(0), fIsVtxErrGuess(0), fIsLinearized(0), fConstructMethod(2), SumDaughterMass(0), fMassHypo(-1), fId(-1), fDaughterIds(), fPDG(0)
28{
29 //* Constructor
30
31 Initialize();
32}
33
34void KFParticleBaseSIMD::Initialize( const fvec Param[], const fvec Cov[], Int_t Charge, fvec 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 fvec energy = 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 fvec energyInv = 1./energy;
63 fvec
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.;
98 fIsVtxGuess = 0;
99 fIsVtxGuess = 0;
100 fIsLinearized = 0;
101 SumDaughterMass = 0;
102 fMassHypo = -1;
103}
104
106{
107 //* Set decay vertex parameters for linearisation
108
109 fVtxGuess[0] = x;
110 fVtxGuess[1] = y;
111 fVtxGuess[2] = z;
112 fIsLinearized = 1;
113}
114
116{
117 //* Set errors of the decay vertex parameters for linearisation
118
119 fVtxErrGuess[0] = dx;
120 fVtxErrGuess[1] = dy;
121 fVtxErrGuess[2] = dz;
122 fIsVtxErrGuess = 1;
123}
124
126{
127 //* Calculate particle momentum
128
129 fvec ret = Zero;
130
131 fvec x = fP[3];
132 fvec y = fP[4];
133 fvec z = fP[5];
134 fvec x2 = x*x;
135 fvec y2 = y*y;
136 fvec z2 = z*z;
137 fvec p2 = x2+y2+z2;
138 p = sqrt(p2);
139 error = (x2*fC[9]+y2*fC[14]+z2*fC[20] + 2*(x*y*fC[13]+x*z*fC[18]+y*z*fC[19]) );
140 const fvec LocalSmall = 1.e-4;
141 fvec mask = fvec( fvec(Zero < error) & fvec(LocalSmall < fabs(p)));
142 error = (mask & error);
143 error = sqrt(error);
144 ret = (!mask);
145 return ret;
146}
147
149{
150 //* Calculate particle transverse momentum
151
152 fvec ret = Zero;
153
154 fvec px = fP[3];
155 fvec py = fP[4];
156 fvec px2 = px*px;
157 fvec py2 = py*py;
158 fvec pt2 = px2+py2;
159 pt = sqrt(pt2);
160 error = (px2*fC[9] + py2*fC[14] + 2*px*py*fC[13] );
161 const fvec LocalSmall = 1.e-8;
162 fvec mask = fvec( fvec(Zero < error) & fvec(LocalSmall < fabs(pt)));
163 error = (mask & error);
164 error = sqrt(error);
165 error = error/pt;
166 error += ((!mask) & fvec(1.e10));
167 ret = (!mask);
168 return ret;
169}
170
172{
173 //* Calculate particle pseudorapidity
174
175 fvec ret = Zero;
176
177 const fvec BIG = 1.e10;
178 const fvec LocalSmall = 1.e-8;
179
180 fvec px = fP[3];
181 fvec py = fP[4];
182 fvec pz = fP[5];
183 fvec pt2 = px*px + py*py;
184 fvec p2 = pt2 + pz*pz;
185 fvec p = sqrt(p2);
186 fvec a = p + pz;
187 fvec b = p - pz;
188 eta = BIG;
189 fvec c = Zero;
190 c = fvec(b > LocalSmall) & (a/b);
191 fvec logc = 0.5*log(c);
192 eta = if3(fvec(LocalSmall<fabs(c)), logc, eta);
193
194 fvec h3 = -px*pz;
195 fvec h4 = -py*pz;
196 fvec pt4 = pt2*pt2;
197 fvec p2pt4 = p2*pt4;
198 error = (h3*h3*fC[9] + h4*h4*fC[14] + pt4*fC[20] + 2*( h3*(h4*fC[13] + fC[18]*pt2) + pt2*h4*fC[19] ) );
199
200 fvec mask = fvec(fvec(LocalSmall < fabs(p2pt4)) & fvec(Zero < error));
201 error = if3( mask, (sqrt(error/p2pt4)), BIG);
202
203 ret = (!mask);
204
205 return ret;
206}
207
209{
210 //* Calculate particle polar angle
211
212 fvec ret = Zero;
213
214 fvec px = fP[3];
215 fvec py = fP[4];
216 fvec px2 = px*px;
217 fvec py2 = py*py;
218 fvec pt2 = px2 + py2;
219 phi = atan2(py,px);
220 error = (py2*fC[9] + px2*fC[14] - fvec(2.f)*px*py*fC[13] );
221
222 fvec mask = fvec(fvec(Zero < error) & fvec(fvec(1.e-4) < pt2));
223 error = if3( mask, sqrt(error)/pt2, fvec(1.e10));
224 ret = !mask;
225 return ret;
226}
227
229{
230 //* Calculate distance to the origin
231
232 fvec ret = Zero;
233
234 fvec x = fP[0];
235 fvec y = fP[1];
236 fvec x2 = x*x;
237 fvec y2 = y*y;
238 r = sqrt(x2 + y2);
239 error = (x2*fC[0] + y2*fC[2] - fvec(2.f)*x*y*fC[1] );
240
241 fvec mask = fvec(fvec(Zero < error) & fvec(fvec(1.e-4) < r));
242 error = if3( mask, sqrt(error)/r, fvec(1.e10));
243 ret = !mask;
244 return ret;
245}
246
248{
249 //* Calculate particle mass
250
251 // s = sigma^2 of m2/2
252
253 const fvec BIG = 1.e20;
254 const fvec LocalSmall = 1.e-10;
255
256 fvec ret = Zero;
257
258 fvec s = ( fP[3]*fP[3]*fC[9] + fP[4]*fP[4]*fC[14] + fP[5]*fP[5]*fC[20]
259 + fP[6]*fP[6]*fC[27]
260 +fvec(2.f)*( + fP[3]*fP[4]*fC[13] + fP[5]*(fP[3]*fC[18] + fP[4]*fC[19])
261 - fP[6]*( fP[3]*fC[24] + fP[4]*fC[25] + fP[5]*fC[26] ) )
262 );
263// fvec m2 = fabs(fP[6]*fP[6] - fP[3]*fP[3] - fP[4]*fP[4] - fP[5]*fP[5]);
264// m = sqrt(m2);
265// if( m>1.e-10 ){
266// if( s>=0 ){
267// error = sqrt(s)/m;
268// return 0;
269// }
270// }
271// error = 1.e20;
272 fvec m2 = (fP[6]*fP[6] - fP[3]*fP[3] - fP[4]*fP[4] - fP[5]*fP[5]);
273
274 fvec mask = fvec(Zero <= m2);
275 m = if3(mask, sqrt(m2), -sqrt(-m2));
276
277 mask = fvec(mask & fvec(fvec(Zero <= s) & (LocalSmall < m)));
278 error = if3(mask, sqrt(s)/m, BIG);
279
280 ret = !mask;
281 return ret;
282}
283
284
286{
287 //* Calculate particle decay length [cm]
288
289 fvec ret = Zero;
290 const fvec BIG = 1.e20;
291
292 fvec x = fP[3];
293 fvec y = fP[4];
294 fvec z = fP[5];
295 fvec t = fP[7];
296 fvec x2 = x*x;
297 fvec y2 = y*y;
298 fvec z2 = z*z;
299 fvec p2 = x2+y2+z2;
300 l = t*sqrt(p2);
301
302 error = p2*fC[35] + t*t/p2*(x2*fC[9]+y2*fC[14]+z2*fC[20]
303 + fvec(2.f)*(x*y*fC[13]+x*z*fC[18]+y*z*fC[19]) )
304 + fvec(2.f)*t*(x*fC[31]+y*fC[32]+z*fC[33]);
305
306 fvec mask = fvec(fvec(1.e-4) < p2);
307 error = if3( mask, sqrt(fabs(error)), BIG);
308 ret = !mask;
309 return ret;
310}
311
313{
314 //* Calculate particle decay length in XY projection [cm]
315
316 const fvec BIG = 1.e20;
317 fvec ret = Zero;
318
319 fvec x = fP[3];
320 fvec y = fP[4];
321 fvec t = fP[7];
322 fvec x2 = x*x;
323 fvec y2 = y*y;
324 fvec pt2 = x2+y2;
325 l = t*sqrt(pt2);
326
327 error = pt2*fC[35] + t*t/pt2*(x2*fC[9]+y2*fC[14] + 2*x*y*fC[13] )
328 + fvec(2.f)*t*(x*fC[31]+y*fC[32]);
329 fvec mask = fvec(fvec(1.e-4) < pt2);
330 error = if3( mask, sqrt(fabs(error)), BIG);
331 ret = !mask;
332 return ret;
333}
334
335
337{
338 //* Calculate particle decay time [s]
339
340 const fvec BIG = 1.e20;
341 fvec ret = Zero;
342
343 fvec m, dm;
344 GetMass( m, dm );
345 fvec cTM = (-fP[3]*fC[31] - fP[4]*fC[32] - fP[5]*fC[33] + fP[6]*fC[34]);
346 tauC = fP[7]*m;
347 error = m*m*fC[35] + 2*fP[7]*cTM + fP[7]*fP[7]*dm*dm;
348 fvec mask = fvec(Zero < error);
349 error = if3( mask, sqrt(error), BIG);
350 ret = !mask;
351 return ret;
352}
353
354
356{
357 //* Add daughter via operator+=
358
359 AddDaughter( Daughter );
360}
361
362fvec KFParticleBaseSIMD::GetSCorrection( const fvec Part[], const fvec XYZ[] ) const
363{
364 //* Get big enough correction for S error to let the particle Part be fitted to XYZ point
365
366 fvec d[3] = { XYZ[0]-Part[0], XYZ[1]-Part[1], XYZ[2]-Part[2] };
367 fvec p2 = Part[3]*Part[3]+Part[4]*Part[4]+Part[5]*Part[5];
368// fvec sigmaS = if3( fvec(fvec(1.e-4)<p2) , ( fvec(10.1)+fvec(3.)*sqrt( d[0]*d[0]+d[1]*d[1]+d[2]*d[2]) )/sqrt(p2) , One);
369// fvec sigmaS = if3( fvec(fvec(1.e-4)<p2) , ( fvec(0.1)+fvec(10.)*sqrt( d[0]*d[0]+d[1]*d[1]+d[2]*d[2]) )/sqrt(p2) , One);
370
371 fvec sigmaS = if3( fvec(fvec(1.e-4)<p2) , fvec(0.1)+fvec(10.)*sqrt( (d[0]*d[0]+d[1]*d[1]+d[2]*d[2])/p2 ) , One);
372// fvec sigmaS = GetDStoPoint(XYZ)*0.08;
373
374 return sigmaS;
375}
376
377void KFParticleBaseSIMD::GetMeasurement( const fvec XYZ[], fvec m[], fvec V[], Bool_t isAtVtxGuess ) const
378{
379 //* Get additional covariances V used during measurement
380
381 fvec b[3];
382 GetFieldValue( XYZ, b );
383 const fvec kCLight = 0.000299792458;
384 b[0]*=kCLight*GetQ(); b[1]*=kCLight*GetQ(); b[2]*=kCLight*GetQ();
385
386 if(!isAtVtxGuess)
387 Transport( GetDStoPoint(XYZ), m, V );
388 else
389 {
390 for(int iM=0; iM<8; iM++)
391 m[iM] = fP[iM];
392 for(int iV=0; iV<8; iV++)
393 V[iV] = fC[iV];
394 }
395
396 fvec sigmaS = GetSCorrection( m, XYZ );
397
398 fvec h[6];
399
400 h[0] = m[3]*sigmaS;
401 h[1] = m[4]*sigmaS;
402 h[2] = m[5]*sigmaS;
403 h[3] = ( h[1]*b[2]-h[2]*b[1] );
404 h[4] = ( h[2]*b[0]-h[0]*b[2] );
405 h[5] = ( h[0]*b[1]-h[1]*b[0] );
406
407 V[ 0]+= h[0]*h[0];
408 V[ 1]+= h[1]*h[0];
409 V[ 2]+= h[1]*h[1];
410 V[ 3]+= h[2]*h[0];
411 V[ 4]+= h[2]*h[1];
412 V[ 5]+= h[2]*h[2];
413
414 V[ 6]+= h[3]*h[0];
415 V[ 7]+= h[3]*h[1];
416 V[ 8]+= h[3]*h[2];
417 V[ 9]+= h[3]*h[3];
418
419 V[10]+= h[4]*h[0];
420 V[11]+= h[4]*h[1];
421 V[12]+= h[4]*h[2];
422 V[13]+= h[4]*h[3];
423 V[14]+= h[4]*h[4];
424
425 V[15]+= h[5]*h[0];
426 V[16]+= h[5]*h[1];
427 V[17]+= h[5]*h[2];
428 V[18]+= h[5]*h[3];
429 V[19]+= h[5]*h[4];
430 V[20]+= h[5]*h[5];
431
432}
433
434void KFParticleBaseSIMD::AddDaughter( const KFParticleBaseSIMD &Daughter, Bool_t isAtVtxGuess )
435{
436 AddDaughterId( Daughter.Id() );
437
438 if( fNDF[0]<-1 ){ // first daughter -> just copy
439 fNDF = -1;
440 fQ = Daughter.GetQ();
441 if( Daughter.fC[35][0]>0 ){ //TODO Check this: only the first daughter is used here!
442 Daughter.GetMeasurement( fVtxGuess, fP, fC, isAtVtxGuess );
443 } else {
444 for( Int_t i=0; i<8; i++ ) fP[i] = Daughter.fP[i];
445 for( Int_t i=0; i<36; i++ ) fC[i] = Daughter.fC[i];
446 }
447 fSFromDecay = 0;
448 fMassHypo = Daughter.fMassHypo;
450 return;
451 }
452
453 if(fConstructMethod == 0)
454 AddDaughterWithEnergyFit(Daughter,isAtVtxGuess);
455 else if(fConstructMethod == 1)
456 AddDaughterWithEnergyCalc(Daughter,isAtVtxGuess);
457 else if(fConstructMethod == 2)
458 AddDaughterWithEnergyFitMC(Daughter,isAtVtxGuess);
459
461 fMassHypo = -1;
462}
463
464void KFParticleBaseSIMD::AddDaughterWithEnergyFit( const KFParticleBaseSIMD &Daughter, Bool_t isAtVtxGuess )
465{
466 //* Energy considered as an independent veriable, fitted independently from momentum, without any constraints on mass
467
468 //* Add daughter
469
470// if(!isAtVtxGuess)
471// TransportToDecayVertex();
472
473 Int_t maxIter = 1;
474
475 if( (!fIsLinearized) && (!isAtVtxGuess) ){
476 if( fNDF[0]==-1 ){
477 fvec ds, ds1;
478 GetDStoParticle(Daughter, ds, ds1);
479 TransportToDS( ds );
480 fvec m[8];
481 fvec mCd[36];
482 Daughter.Transport( ds1, m, mCd );
483 fVtxGuess[0] = .5*( fP[0] + m[0] );
484 fVtxGuess[1] = .5*( fP[1] + m[1] );
485 fVtxGuess[2] = .5*( fP[2] + m[2] );
486 } else {
487 fVtxGuess[0] = fP[0];
488 fVtxGuess[1] = fP[1];
489 fVtxGuess[2] = fP[2];
490 }
491 maxIter = 3;
492 }
493
494 for( Int_t iter=0; iter<maxIter; iter++ ){
495
496 fvec *ffP = fP, *ffC = fC;
497
498 fvec m[8], mV[36];
499
500 if( Daughter.fC[35][0]>0 ){ //TODO Check this: only the first daughter is used here!
501 Daughter.GetMeasurement( fVtxGuess, m, mV, isAtVtxGuess );
502 } else {
503 for( Int_t i=0; i<8; i++ ) m[i] = Daughter.fP[i];
504 for( Int_t i=0; i<36; i++ ) mV[i] = Daughter.fC[i];
505 }
506 //*
507
508 fvec mS[6]= { ffC[0]+mV[0],
509 ffC[1]+mV[1], ffC[2]+mV[2],
510 ffC[3]+mV[3], ffC[4]+mV[4], ffC[5]+mV[5] };
512
513 //* Residual (measured - estimated)
514
515 fvec zeta[3] = { m[0]-ffP[0], m[1]-ffP[1], m[2]-ffP[2] };
516
517
518 //* CHt = CH' - D'
519
520 fvec mCHt0[7], mCHt1[7], mCHt2[7];
521
522 mCHt0[0]=ffC[ 0] ; mCHt1[0]=ffC[ 1] ; mCHt2[0]=ffC[ 3] ;
523 mCHt0[1]=ffC[ 1] ; mCHt1[1]=ffC[ 2] ; mCHt2[1]=ffC[ 4] ;
524 mCHt0[2]=ffC[ 3] ; mCHt1[2]=ffC[ 4] ; mCHt2[2]=ffC[ 5] ;
525 mCHt0[3]=ffC[ 6]-mV[ 6]; mCHt1[3]=ffC[ 7]-mV[ 7]; mCHt2[3]=ffC[ 8]-mV[ 8];
526 mCHt0[4]=ffC[10]-mV[10]; mCHt1[4]=ffC[11]-mV[11]; mCHt2[4]=ffC[12]-mV[12];
527 mCHt0[5]=ffC[15]-mV[15]; mCHt1[5]=ffC[16]-mV[16]; mCHt2[5]=ffC[17]-mV[17];
528 mCHt0[6]=ffC[21]-mV[21]; mCHt1[6]=ffC[22]-mV[22]; mCHt2[6]=ffC[23]-mV[23];
529
530 //* Kalman gain K = mCH'*S
531
532 fvec k0[7], k1[7], k2[7];
533
534 for(Int_t i=0;i<7;++i){
535 k0[i] = mCHt0[i]*mS[0] + mCHt1[i]*mS[1] + mCHt2[i]*mS[3];
536 k1[i] = mCHt0[i]*mS[1] + mCHt1[i]*mS[2] + mCHt2[i]*mS[4];
537 k2[i] = mCHt0[i]*mS[3] + mCHt1[i]*mS[4] + mCHt2[i]*mS[5];
538 }
539
540 //* New estimation of the vertex position
541
542 if( iter<maxIter-1 ){
543 for(Int_t i=0; i<3; ++i)
544 fVtxGuess[i]= ffP[i] + k0[i]*zeta[0]+k1[i]*zeta[1]+k2[i]*zeta[2];
545 continue;
546 }
547
548 // last itearation -> update the particle
549
550 //* Add the daughter momentum to the particle momentum
551
552 ffP[ 3] += m[ 3];
553 ffP[ 4] += m[ 4];
554 ffP[ 5] += m[ 5];
555 ffP[ 6] += m[ 6];
556
557 ffC[ 9] += mV[ 9];
558 ffC[13] += mV[13];
559 ffC[14] += mV[14];
560 ffC[18] += mV[18];
561 ffC[19] += mV[19];
562 ffC[20] += mV[20];
563 ffC[24] += mV[24];
564 ffC[25] += mV[25];
565 ffC[26] += mV[26];
566 ffC[27] += mV[27];
567
568
569 //* New estimation of the vertex position r += K*zeta
570
571 for(Int_t i=0;i<7;++i)
572 fP[i] = ffP[i] + (k0[i]*zeta[0] + k1[i]*zeta[1] + k2[i]*zeta[2]);
573
574 //* New covariance matrix C -= K*(mCH')'
575
576 for(Int_t i=0, k=0;i<7;++i){
577 for(Int_t j=0;j<=i;++j,++k){
578 fC[k] = ffC[k] - (k0[i]*mCHt0[j] + k1[i]*mCHt1[j] + k2[i]*mCHt2[j] );
579 }
580 }
581
582 //* Calculate Chi^2
583 if( iter == (maxIter-1) ) {
584 fNDF += 2;
585 fQ += Daughter.GetQ();
586 fSFromDecay = 0;
587 fChi2 += (mS[0]*zeta[0] + mS[1]*zeta[1] + mS[3]*zeta[2])*zeta[0]
588 + (mS[1]*zeta[0] + mS[2]*zeta[1] + mS[4]*zeta[2])*zeta[1]
589 + (mS[3]*zeta[0] + mS[4]*zeta[1] + mS[5]*zeta[2])*zeta[2];
590 }
591
592 }
593}
594
595void KFParticleBaseSIMD::AddDaughterWithEnergyCalc( const KFParticleBaseSIMD &Daughter, Bool_t isAtVtxGuess )
596{
597 //* Energy considered as a dependent variable, calculated from the momentum and mass hypothesis
598
599 //* Add daughter
600
601// if(!isAtVtxGuess)
602// TransportToDecayVertex();
603
604 Int_t maxIter = 1;
605
606 if( (!fIsLinearized) && (!isAtVtxGuess) ){
607 if( fNDF[0]==-1 ){
608 fvec ds, ds1;
609 GetDStoParticle(Daughter, ds, ds1);
610 TransportToDS( ds );
611 fvec m[8];
612 fvec mCd[36];
613 Daughter.Transport( ds1, m, mCd );
614 fVtxGuess[0] = .5*( fP[0] + m[0] );
615 fVtxGuess[1] = .5*( fP[1] + m[1] );
616 fVtxGuess[2] = .5*( fP[2] + m[2] );
617 } else {
618 fVtxGuess[0] = fP[0];
619 fVtxGuess[1] = fP[1];
620 fVtxGuess[2] = fP[2];
621 }
622 maxIter = 3;
623 }
624
625 for( Int_t iter=0; iter<maxIter; iter++ ){
626
627 fvec *ffP = fP, *ffC = fC;
628
629 fvec m[8], mV[36];
630
631 if( Daughter.fC[35][0]>0 ){ //TODO Check this: only the first daughter is used here!
632 Daughter.GetMeasurement( fVtxGuess, m, mV, isAtVtxGuess );
633 } else {
634 for( Int_t i=0; i<8; i++ ) m[i] = Daughter.fP[i];
635 for( Int_t i=0; i<36; i++ ) mV[i] = Daughter.fC[i];
636 }
637
638 fvec massMf2 = m[6]*m[6] - (m[3]*m[3] + m[4]*m[4] + m[5]*m[5]);
639 fvec massRf2 = fP[6]*fP[6] - (fP[3]*fP[3] + fP[4]*fP[4] + fP[5]*fP[5]);
640
641 //*
642
643 fvec mS[6]= { ffC[0]+mV[0],
644 ffC[1]+mV[1], ffC[2]+mV[2],
645 ffC[3]+mV[3], ffC[4]+mV[4], ffC[5]+mV[5] };
647
648 //* Residual (measured - estimated)
649
650 fvec zeta[3] = { m[0]-ffP[0], m[1]-ffP[1], m[2]-ffP[2] };
651
652 //* CHt = CH' - D'
653
654 fvec mCHt0[6], mCHt1[6], mCHt2[6];
655
656 mCHt0[0]=ffC[ 0] ; mCHt1[0]=ffC[ 1] ; mCHt2[0]=ffC[ 3] ;
657 mCHt0[1]=ffC[ 1] ; mCHt1[1]=ffC[ 2] ; mCHt2[1]=ffC[ 4] ;
658 mCHt0[2]=ffC[ 3] ; mCHt1[2]=ffC[ 4] ; mCHt2[2]=ffC[ 5] ;
659 mCHt0[3]=ffC[ 6]-mV[ 6]; mCHt1[3]=ffC[ 7]-mV[ 7]; mCHt2[3]=ffC[ 8]-mV[ 8];
660 mCHt0[4]=ffC[10]-mV[10]; mCHt1[4]=ffC[11]-mV[11]; mCHt2[4]=ffC[12]-mV[12];
661 mCHt0[5]=ffC[15]-mV[15]; mCHt1[5]=ffC[16]-mV[16]; mCHt2[5]=ffC[17]-mV[17];
662
663 //* Kalman gain K = mCH'*S
664
665 fvec k0[6], k1[6], k2[6];
666
667 for(Int_t i=0;i<6;++i){
668 k0[i] = mCHt0[i]*mS[0] + mCHt1[i]*mS[1] + mCHt2[i]*mS[3];
669 k1[i] = mCHt0[i]*mS[1] + mCHt1[i]*mS[2] + mCHt2[i]*mS[4];
670 k2[i] = mCHt0[i]*mS[3] + mCHt1[i]*mS[4] + mCHt2[i]*mS[5];
671 }
672
673 //* New estimation of the vertex position
674
675 if( iter<maxIter-1 ){
676 for(Int_t i=0; i<3; ++i)
677 fVtxGuess[i]= ffP[i] + k0[i]*zeta[0]+k1[i]*zeta[1]+k2[i]*zeta[2];
678 continue;
679 }
680
681 //* find mf and mVf - optimum value of the measurement and its covariance matrix
682 //* mVHt = V*H'
683 fvec mVHt0[6], mVHt1[6], mVHt2[6];
684
685 mVHt0[0]= mV[ 0] ; mVHt1[0]= mV[ 1] ; mVHt2[0]= mV[ 3] ;
686 mVHt0[1]= mV[ 1] ; mVHt1[1]= mV[ 2] ; mVHt2[1]= mV[ 4] ;
687 mVHt0[2]= mV[ 3] ; mVHt1[2]= mV[ 4] ; mVHt2[2]= mV[ 5] ;
688 mVHt0[3]= mV[ 6] ; mVHt1[3]= mV[ 7] ; mVHt2[3]= mV[ 8] ;
689 mVHt0[4]= mV[10] ; mVHt1[4]= mV[11] ; mVHt2[4]= mV[12] ;
690 mVHt0[5]= mV[15] ; mVHt1[5]= mV[16] ; mVHt2[5]= mV[17] ;
691
692 //* Kalman gain Km = mCH'*S
693
694 fvec km0[6], km1[6], km2[6];
695
696 for(Int_t i=0;i<6;++i){
697 km0[i] = mVHt0[i]*mS[0] + mVHt1[i]*mS[1] + mVHt2[i]*mS[3];
698 km1[i] = mVHt0[i]*mS[1] + mVHt1[i]*mS[2] + mVHt2[i]*mS[4];
699 km2[i] = mVHt0[i]*mS[3] + mVHt1[i]*mS[4] + mVHt2[i]*mS[5];
700 }
701
702 fvec mf[7] = { m[0], m[1], m[2], m[3], m[4], m[5], m[6] };
703
704 for(Int_t i=0;i<6;++i)
705 mf[i] = mf[i] - km0[i]*zeta[0] - km1[i]*zeta[1] - km2[i]*zeta[2];
706
707 fvec energyMf = sqrt( massMf2 + (mf[3]*mf[3] + mf[4]*mf[4] + mf[5]*mf[5]) );
708
709 fvec mVf[28];
710 for(Int_t iC=0; iC<28; iC++)
711 mVf[iC] = mV[iC];
712
713 //* hmf = d(energyMf)/d(mf)
714 fvec hmf[7];
715 hmf[3] = if3( fvec(fabs(energyMf) < small), Zero, hmf[3]/energyMf);
716 hmf[4] = if3( fvec(fabs(energyMf) < small), Zero, hmf[4]/energyMf);
717 hmf[5] = if3( fvec(fabs(energyMf) < small), Zero, hmf[5]/energyMf);
718 hmf[6] = 0;
719
720 for(Int_t i=0, k=0;i<6;++i){
721 for(Int_t j=0;j<=i;++j,++k){
722 mVf[k] = mVf[k] - (km0[i]*mVHt0[j] + km1[i]*mVHt1[j] + km2[i]*mVHt2[j] );
723 }
724 }
725 fvec mVf24 = mVf[24], mVf25 = mVf[25], mVf26 = mVf[26];
726 mVf[21] = mVf[6 ]*hmf[3] + mVf[10]*hmf[4] + mVf[15]*hmf[5] + mVf[21]*hmf[6];
727 mVf[22] = mVf[7 ]*hmf[3] + mVf[11]*hmf[4] + mVf[16]*hmf[5] + mVf[22]*hmf[6];
728 mVf[23] = mVf[8 ]*hmf[3] + mVf[12]*hmf[4] + mVf[17]*hmf[5] + mVf[23]*hmf[6];
729 mVf[24] = mVf[9 ]*hmf[3] + mVf[13]*hmf[4] + mVf[18]*hmf[5] + mVf[24]*hmf[6];
730 mVf[25] = mVf[13]*hmf[3] + mVf[14]*hmf[4] + mVf[19]*hmf[5] + mVf[25]*hmf[6];
731 mVf[26] = mVf[18]*hmf[3] + mVf[19]*hmf[4] + mVf[20]*hmf[5] + mVf[26]*hmf[6];
732 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
733
734 mf[6] = energyMf;
735
736 //* find rf and mCf - optimum value of the measurement and its covariance matrix
737
738 //* mCCHt = C*H'
739 fvec mCCHt0[6], mCCHt1[6], mCCHt2[6];
740
741 mCCHt0[0]=ffC[ 0]; mCCHt1[0]=ffC[ 1]; mCCHt2[0]=ffC[ 3];
742 mCCHt0[1]=ffC[ 1]; mCCHt1[1]=ffC[ 2]; mCCHt2[1]=ffC[ 4];
743 mCCHt0[2]=ffC[ 3]; mCCHt1[2]=ffC[ 4]; mCCHt2[2]=ffC[ 5];
744 mCCHt0[3]=ffC[ 6]; mCCHt1[3]=ffC[ 7]; mCCHt2[3]=ffC[ 8];
745 mCCHt0[4]=ffC[10]; mCCHt1[4]=ffC[11]; mCCHt2[4]=ffC[12];
746 mCCHt0[5]=ffC[15]; mCCHt1[5]=ffC[16]; mCCHt2[5]=ffC[17];
747
748 //* Kalman gain Krf = mCH'*S
749
750 fvec krf0[6], krf1[6], krf2[6];
751
752 for(Int_t i=0;i<6;++i){
753 krf0[i] = mCCHt0[i]*mS[0] + mCCHt1[i]*mS[1] + mCCHt2[i]*mS[3];
754 krf1[i] = mCCHt0[i]*mS[1] + mCCHt1[i]*mS[2] + mCCHt2[i]*mS[4];
755 krf2[i] = mCCHt0[i]*mS[3] + mCCHt1[i]*mS[4] + mCCHt2[i]*mS[5];
756 }
757 fvec rf[7] = { ffP[0], ffP[1], ffP[2], ffP[3], ffP[4], ffP[5], ffP[6] };
758
759 for(Int_t i=0;i<6;++i)
760 rf[i] = rf[i] + krf0[i]*zeta[0] + krf1[i]*zeta[1] + krf2[i]*zeta[2];
761
762 fvec energyRf = sqrt( massRf2 + (rf[3]*rf[3] + rf[4]*rf[4] + rf[5]*rf[5]) );
763
764 fvec mCf[28];
765 for(Int_t iC=0; iC<28; iC++)
766 mCf[iC] = ffC[iC];
767 //* hrf = d(Erf)/d(rf)
768 fvec hrf[7];
769 hrf[3] = if3( fvec(fabs(energyRf) < small), Zero, rf[3]/energyRf);
770 hrf[4] = if3( fvec(fabs(energyRf) < small), Zero, rf[4]/energyRf);
771 hrf[5] = if3( fvec(fabs(energyRf) < small), Zero, rf[5]/energyRf);
772// hrf[6] = if3( fvec(fabs(E_rf) < small), Zero, rf[6]/E_rf);
773 hrf[6] = 0;
774
775 for(Int_t i=0, k=0;i<6;++i){
776 for(Int_t j=0;j<=i;++j,++k){
777 mCf[k] = mCf[k] - (krf0[i]*mCCHt0[j] + krf1[i]*mCCHt1[j] + krf2[i]*mCCHt2[j] );
778 }
779 }
780 fvec mCf24 = mCf[24], mCf25 = mCf[25], mCf26 = mCf[26];
781 mCf[21] = mCf[6 ]*hrf[3] + mCf[10]*hrf[4] + mCf[15]*hrf[5] + mCf[21]*hrf[6];
782 mCf[22] = mCf[7 ]*hrf[3] + mCf[11]*hrf[4] + mCf[16]*hrf[5] + mCf[22]*hrf[6];
783 mCf[23] = mCf[8 ]*hrf[3] + mCf[12]*hrf[4] + mCf[17]*hrf[5] + mCf[23]*hrf[6];
784 mCf[24] = mCf[9 ]*hrf[3] + mCf[13]*hrf[4] + mCf[18]*hrf[5] + mCf[24]*hrf[6];
785 mCf[25] = mCf[13]*hrf[3] + mCf[14]*hrf[4] + mCf[19]*hrf[5] + mCf[25]*hrf[6];
786 mCf[26] = mCf[18]*hrf[3] + mCf[19]*hrf[4] + mCf[20]*hrf[5] + mCf[26]*hrf[6];
787 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
788
789 for(Int_t iC=21; iC<28; iC++)
790 {
791 ffC[iC] = mCf[iC];
792 mV[iC] = mVf[iC];
793 }
794
795 fP[6] = energyRf + energyMf;
796 rf[6] = energyRf;
797
798 //fvec Dvv[3][3]; do not need this
799 fvec mDvp[3][3];
800 fvec mDpv[3][3];
801 fvec mDpp[3][3];
802 fvec mDe[7];
803
804 for(int i=0; i<3; i++)
805 {
806 for(int j=0; j<3; j++)
807 {
808 mDvp[i][j] = km0[i+3]*mCCHt0[j] + km1[i+3]*mCCHt1[j] + km2[i+3]*mCCHt2[j];
809 mDpv[i][j] = km0[i]*mCCHt0[j+3] + km1[i]*mCCHt1[j+3] + km2[i]*mCCHt2[j+3];
810 mDpp[i][j] = km0[i+3]*mCCHt0[j+3] + km1[i+3]*mCCHt1[j+3] + km2[i+3]*mCCHt2[j+3];
811 }
812 }
813
814 mDe[0] = hmf[3]*mDvp[0][0] + hmf[4]*mDvp[1][0] + hmf[5]*mDvp[2][0];
815 mDe[1] = hmf[3]*mDvp[0][1] + hmf[4]*mDvp[1][1] + hmf[5]*mDvp[2][1];
816 mDe[2] = hmf[3]*mDvp[0][2] + hmf[4]*mDvp[1][2] + hmf[5]*mDvp[2][2];
817 mDe[3] = hmf[3]*mDpp[0][0] + hmf[4]*mDpp[1][0] + hmf[5]*mDpp[2][0];
818 mDe[4] = hmf[3]*mDpp[0][1] + hmf[4]*mDpp[1][1] + hmf[5]*mDpp[2][1];
819 mDe[5] = hmf[3]*mDpp[0][2] + hmf[4]*mDpp[1][2] + hmf[5]*mDpp[2][2];
820 mDe[6] = 2*(mDe[3]*hrf[3] + mDe[4]*hrf[4] + mDe[5]*hrf[5]);
821
822 // last itearation -> update the particle
823
824 //* Add the daughter momentum to the particle momentum
825
826 ffP[ 3] += m[ 3];
827 ffP[ 4] += m[ 4];
828 ffP[ 5] += m[ 5];
829
830 ffC[ 9] += mV[ 9];
831 ffC[13] += mV[13];
832 ffC[14] += mV[14];
833 ffC[18] += mV[18];
834 ffC[19] += mV[19];
835 ffC[20] += mV[20];
836 ffC[24] += mV[24];
837 ffC[25] += mV[25];
838 ffC[26] += mV[26];
839 ffC[27] += mV[27];
840
841 ffC[21] += mDe[0];
842 ffC[22] += mDe[1];
843 ffC[23] += mDe[2];
844 ffC[24] += mDe[3];
845 ffC[25] += mDe[4];
846 ffC[26] += mDe[5];
847 ffC[27] += mDe[6];
848
849 //* New estimation of the vertex position r += K*zeta
850
851 for(Int_t i=0;i<6;++i)
852 fP[i] = ffP[i] + k0[i]*zeta[0] + k1[i]*zeta[1] + k2[i]*zeta[2];
853
854 //* New covariance matrix C -= K*(mCH')'
855
856 for(Int_t i=0, k=0;i<6;++i){
857 for(Int_t j=0;j<=i;++j,++k){
858 fC[k] = ffC[k] - (k0[i]*mCHt0[j] + k1[i]*mCHt1[j] + k2[i]*mCHt2[j] );
859 }
860 }
861
862 for(int i=21; i<28; i++) fC[i] = ffC[i];
863
864 //* Calculate Chi^2
865
866 if( iter == (maxIter-1) ) {
867 fNDF += 2;
868 fQ += Daughter.GetQ();
869 fSFromDecay = 0;
870 fChi2 += (mS[0]*zeta[0] + mS[1]*zeta[1] + mS[3]*zeta[2])*zeta[0]
871 + (mS[1]*zeta[0] + mS[2]*zeta[1] + mS[4]*zeta[2])*zeta[1]
872 + (mS[3]*zeta[0] + mS[4]*zeta[1] + mS[5]*zeta[2])*zeta[2];
873 }
874 }
875}
876
877void KFParticleBaseSIMD::AddDaughterWithEnergyFitMC( const KFParticleBaseSIMD &Daughter, Bool_t isAtVtxGuess )
878{
879 //* Energy considered as an independent variable, fitted independently from momentum, without any constraints on mass
880
881 //* Add daughter
882
883// if(!isAtVtxGuess)
884// TransportToDecayVertex();
885
886 Int_t maxIter = 1;
887
888 if( (!fIsLinearized) && (!isAtVtxGuess) ){
889 if( fNDF[0]==-1 ){
890 fvec ds, ds1;
891 GetDStoParticle(Daughter, ds, ds1);
892 TransportToDS( ds );
893 fvec m[8];
894 fvec mCd[36];
895 Daughter.Transport( ds1, m, mCd );
896 fVtxGuess[0] = .5*( fP[0] + m[0] );
897 fVtxGuess[1] = .5*( fP[1] + m[1] );
898 fVtxGuess[2] = .5*( fP[2] + m[2] );
899 } else {
900 fVtxGuess[0] = fP[0];
901 fVtxGuess[1] = fP[1];
902 fVtxGuess[2] = fP[2];
903 }
904 maxIter = 3;
905 }
906
907 for( Int_t iter=0; iter<maxIter; iter++ ){
908
909 fvec *ffP = fP, *ffC = fC;
910 fvec m[8], mV[36];
911
912 if( Daughter.fC[35][0]>0 ){ //TODO Check this: only the first daughter is used here!
913 Daughter.GetMeasurement( fVtxGuess, m, mV, isAtVtxGuess );
914 } else {
915 for( Int_t i=0; i<8; i++ ) m[i] = Daughter.fP[i];
916 for( Int_t i=0; i<36; i++ ) mV[i] = Daughter.fC[i];
917 }
918 //*
919
920 fvec mS[6]= { ffC[0]+mV[0],
921 ffC[1]+mV[1], ffC[2]+mV[2],
922 ffC[3]+mV[3], ffC[4]+mV[4], ffC[5]+mV[5] };
924
925 //* Residual (measured - estimated)
926
927 fvec zeta[3] = { m[0]-ffP[0], m[1]-ffP[1], m[2]-ffP[2] };
928
929
930 //* CHt = CH'
931
932 fvec mCHt0[7], mCHt1[7], mCHt2[7];
933
934 mCHt0[0]=ffC[ 0] ; mCHt1[0]=ffC[ 1] ; mCHt2[0]=ffC[ 3] ;
935 mCHt0[1]=ffC[ 1] ; mCHt1[1]=ffC[ 2] ; mCHt2[1]=ffC[ 4] ;
936 mCHt0[2]=ffC[ 3] ; mCHt1[2]=ffC[ 4] ; mCHt2[2]=ffC[ 5] ;
937 mCHt0[3]=ffC[ 6] ; mCHt1[3]=ffC[ 7] ; mCHt2[3]=ffC[ 8] ;
938 mCHt0[4]=ffC[10] ; mCHt1[4]=ffC[11] ; mCHt2[4]=ffC[12] ;
939 mCHt0[5]=ffC[15] ; mCHt1[5]=ffC[16] ; mCHt2[5]=ffC[17] ;
940 mCHt0[6]=ffC[21] ; mCHt1[6]=ffC[22] ; mCHt2[6]=ffC[23] ;
941
942 //* Kalman gain K = mCH'*S
943
944 fvec k0[7], k1[7], k2[7];
945
946 for(Int_t i=0;i<7;++i){
947 k0[i] = mCHt0[i]*mS[0] + mCHt1[i]*mS[1] + mCHt2[i]*mS[3];
948 k1[i] = mCHt0[i]*mS[1] + mCHt1[i]*mS[2] + mCHt2[i]*mS[4];
949 k2[i] = mCHt0[i]*mS[3] + mCHt1[i]*mS[4] + mCHt2[i]*mS[5];
950 }
951
952 //* New estimation of the vertex position
953
954 if( iter<maxIter-1 ){
955 for(Int_t i=0; i<3; ++i)
956 fVtxGuess[i]= ffP[i] + k0[i]*zeta[0]+k1[i]*zeta[1]+k2[i]*zeta[2];
957 continue;
958 }
959
960 // last itearation -> update the particle
961
962 //* VHt = VH'
963
964 fvec mVHt0[7], mVHt1[7], mVHt2[7];
965
966 mVHt0[0]=mV[ 0] ; mVHt1[0]=mV[ 1] ; mVHt2[0]=mV[ 3] ;
967 mVHt0[1]=mV[ 1] ; mVHt1[1]=mV[ 2] ; mVHt2[1]=mV[ 4] ;
968 mVHt0[2]=mV[ 3] ; mVHt1[2]=mV[ 4] ; mVHt2[2]=mV[ 5] ;
969 mVHt0[3]=mV[ 6] ; mVHt1[3]=mV[ 7] ; mVHt2[3]=mV[ 8] ;
970 mVHt0[4]=mV[10] ; mVHt1[4]=mV[11] ; mVHt2[4]=mV[12] ;
971 mVHt0[5]=mV[15] ; mVHt1[5]=mV[16] ; mVHt2[5]=mV[17] ;
972 mVHt0[6]=mV[21] ; mVHt1[6]=mV[22] ; mVHt2[6]=mV[23] ;
973
974 //* Kalman gain Km = mCH'*S
975
976 fvec km0[7], km1[7], km2[7];
977
978 for(Int_t i=0;i<7;++i){
979 km0[i] = mVHt0[i]*mS[0] + mVHt1[i]*mS[1] + mVHt2[i]*mS[3];
980 km1[i] = mVHt0[i]*mS[1] + mVHt1[i]*mS[2] + mVHt2[i]*mS[4];
981 km2[i] = mVHt0[i]*mS[3] + mVHt1[i]*mS[4] + mVHt2[i]*mS[5];
982 }
983
984 for(Int_t i=0;i<7;++i)
985 ffP[i] = ffP[i] + k0[i]*zeta[0] + k1[i]*zeta[1] + k2[i]*zeta[2];
986
987 for(Int_t i=0;i<7;++i)
988 m[i] = m[i] - km0[i]*zeta[0] - km1[i]*zeta[1] - km2[i]*zeta[2];
989
990 for(Int_t i=0, k=0;i<7;++i){
991 for(Int_t j=0;j<=i;++j,++k){
992 ffC[k] = ffC[k] - (k0[i]*mCHt0[j] + k1[i]*mCHt1[j] + k2[i]*mCHt2[j] );
993 }
994 }
995
996 for(Int_t i=0, k=0;i<7;++i){
997 for(Int_t j=0;j<=i;++j,++k){
998 mV[k] = mV[k] - (km0[i]*mVHt0[j] + km1[i]*mVHt1[j] + km2[i]*mVHt2[j] );
999 }
1000 }
1001
1002 fvec mDf[7][7];
1003
1004 for(Int_t i=0;i<7;++i){
1005 for(Int_t j=0;j<7;++j){
1006 mDf[i][j] = (km0[i]*mCHt0[j] + km1[i]*mCHt1[j] + km2[i]*mCHt2[j] );
1007 }
1008 }
1009
1010 fvec mJ1[7][7], mJ2[7][7];
1011 for(Int_t iPar1=0; iPar1<7; iPar1++)
1012 {
1013 for(Int_t iPar2=0; iPar2<7; iPar2++)
1014 {
1015 mJ1[iPar1][iPar2] = 0;
1016 mJ2[iPar1][iPar2] = 0;
1017 }
1018 }
1019
1020 fvec mMassParticle = ffP[6]*ffP[6] - (ffP[3]*ffP[3] + ffP[4]*ffP[4] + ffP[5]*ffP[5]);
1021 fvec mMassDaughter = m[6]*m[6] - (m[3]*m[3] + m[4]*m[4] + m[5]*m[5]);
1022 mMassParticle = if3( fvec( mMassParticle > Zero ), sqrt(mMassParticle) , Zero);
1023 mMassDaughter = if3( fvec( mMassDaughter > Zero ), sqrt(mMassDaughter) , Zero);
1024
1025 fvec mask1 = fvec(fMassHypo > fvec(-0.5));
1026 fvec mask2 = fvec ( fvec(!mask1) & fvec( fvec(mMassParticle < SumDaughterMass) | fvec(ffP[6]<Zero)) );
1027 SetMassConstraint(ffP,ffC,mJ1,fMassHypo, mask1);
1028 SetMassConstraint(ffP,ffC,mJ1,SumDaughterMass, mask2);
1029
1030 fvec mask3 = fvec(fvec(Daughter.fMassHypo) > fvec(-0.5));
1031 fvec mask4 = fvec( fvec(!mask3) & fvec( fvec(mMassDaughter<Daughter.SumDaughterMass) | fvec(m[6]<Zero)) );
1032 SetMassConstraint(m,mV,mJ2,Daughter.fMassHypo, mask3);
1033 SetMassConstraint(m,mV,mJ2,Daughter.SumDaughterMass, mask4);
1034
1035 fvec mDJ[7][7];
1036
1037 for(Int_t i=0; i<7; i++) {
1038 for(Int_t j=0; j<7; j++) {
1039 mDJ[i][j] = 0;
1040 for(Int_t k=0; k<7; k++) {
1041 mDJ[i][j] += mDf[i][k]*mJ1[j][k];
1042 }
1043 }
1044 }
1045
1046 for(Int_t i=0; i<7; ++i){
1047 for(Int_t j=0; j<7; ++j){
1048 mDf[i][j]=0;
1049 for(Int_t l=0; l<7; l++){
1050 mDf[i][j] += mJ2[i][l]*mDJ[l][j];
1051 }
1052 }
1053 }
1054
1055 //* Add the daughter momentum to the particle momentum
1056
1057 ffP[ 3] += m[ 3];
1058 ffP[ 4] += m[ 4];
1059 ffP[ 5] += m[ 5];
1060 ffP[ 6] += m[ 6];
1061
1062 ffC[ 9] += mV[ 9];
1063 ffC[13] += mV[13];
1064 ffC[14] += mV[14];
1065 ffC[18] += mV[18];
1066 ffC[19] += mV[19];
1067 ffC[20] += mV[20];
1068 ffC[24] += mV[24];
1069 ffC[25] += mV[25];
1070 ffC[26] += mV[26];
1071 ffC[27] += mV[27];
1072
1073 ffC[6 ] += mDf[3][0]; ffC[7 ] += mDf[3][1]; ffC[8 ] += mDf[3][2];
1074 ffC[10] += mDf[4][0]; ffC[11] += mDf[4][1]; ffC[12] += mDf[4][2];
1075 ffC[15] += mDf[5][0]; ffC[16] += mDf[5][1]; ffC[17] += mDf[5][2];
1076 ffC[21] += mDf[6][0]; ffC[22] += mDf[6][1]; ffC[23] += mDf[6][2];
1077
1078 ffC[9 ] += mDf[3][3] + mDf[3][3];
1079 ffC[13] += mDf[4][3] + mDf[3][4]; ffC[14] += mDf[4][4] + mDf[4][4];
1080 ffC[18] += mDf[5][3] + mDf[3][5]; ffC[19] += mDf[5][4] + mDf[4][5]; ffC[20] += mDf[5][5] + mDf[5][5];
1081 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];
1082
1083 //* New estimation of the vertex position r += K*zeta
1084
1085 for(Int_t i=0;i<7;++i)
1086 fP[i] = ffP[i];
1087
1088 //* New covariance matrix C -= K*(mCH')'
1089
1090 for(Int_t i=0, k=0;i<7;++i){
1091 for(Int_t j=0;j<=i;++j,++k){
1092 fC[k] = ffC[k];
1093 }
1094 }
1095 //* Calculate Chi^2
1096
1097 if( iter == (maxIter-1) ) {
1098 fNDF += 2;
1099 fQ += Daughter.GetQ();
1100 fSFromDecay = 0;
1101 fChi2 += (mS[0]*zeta[0] + mS[1]*zeta[1] + mS[3]*zeta[2])*zeta[0]
1102 + (mS[1]*zeta[0] + mS[2]*zeta[1] + mS[4]*zeta[2])*zeta[1]
1103 + (mS[3]*zeta[0] + mS[4]*zeta[1] + mS[5]*zeta[2])*zeta[2];
1104 }
1105 }
1106}
1107
1109{
1110 //* Set production vertex for the particle, when the particle was not used in the vertex fit
1111
1112 const fvec *m = Vtx.fP, *mV = Vtx.fC;
1113
1114 Bool_t noS = ( fC[35][0]<=0 ); // no decay length allowed //TODO Check this: only the first daughter is used here!
1115
1116 if( noS ){
1118 fP[7] = 0;
1119 fC[28] = fC[29] = fC[30] = fC[31] = fC[32] = fC[33] = fC[34] = fC[35] = 0;
1120 } else {
1122 fP[7] = -fSFromDecay;
1123 fC[28] = fC[29] = fC[30] = fC[31] = fC[32] = fC[33] = fC[34] = 0;
1124 fC[35] = 0.1;
1125
1126 Convert(1);
1127 }
1128
1129 fvec mAi[6] = {fC[0], fC[1], fC[2], fC[3], fC[4], fC[5]};
1130
1131 InvertCholetsky3(mAi);
1132// InvertSym3( fC, mAi );
1133
1134 fvec mB[5][3];
1135
1136 mB[0][0] = fC[ 6]*mAi[0] + fC[ 7]*mAi[1] + fC[ 8]*mAi[3];
1137 mB[0][1] = fC[ 6]*mAi[1] + fC[ 7]*mAi[2] + fC[ 8]*mAi[4];
1138 mB[0][2] = fC[ 6]*mAi[3] + fC[ 7]*mAi[4] + fC[ 8]*mAi[5];
1139
1140 mB[1][0] = fC[10]*mAi[0] + fC[11]*mAi[1] + fC[12]*mAi[3];
1141 mB[1][1] = fC[10]*mAi[1] + fC[11]*mAi[2] + fC[12]*mAi[4];
1142 mB[1][2] = fC[10]*mAi[3] + fC[11]*mAi[4] + fC[12]*mAi[5];
1143
1144 mB[2][0] = fC[15]*mAi[0] + fC[16]*mAi[1] + fC[17]*mAi[3];
1145 mB[2][1] = fC[15]*mAi[1] + fC[16]*mAi[2] + fC[17]*mAi[4];
1146 mB[2][2] = fC[15]*mAi[3] + fC[16]*mAi[4] + fC[17]*mAi[5];
1147
1148 mB[3][0] = fC[21]*mAi[0] + fC[22]*mAi[1] + fC[23]*mAi[3];
1149 mB[3][1] = fC[21]*mAi[1] + fC[22]*mAi[2] + fC[23]*mAi[4];
1150 mB[3][2] = fC[21]*mAi[3] + fC[22]*mAi[4] + fC[23]*mAi[5];
1151
1152 mB[4][0] = fC[28]*mAi[0] + fC[29]*mAi[1] + fC[30]*mAi[3];
1153 mB[4][1] = fC[28]*mAi[1] + fC[29]*mAi[2] + fC[30]*mAi[4];
1154 mB[4][2] = fC[28]*mAi[3] + fC[29]*mAi[4] + fC[30]*mAi[5];
1155
1156 fvec z[3] = { m[0]-fP[0], m[1]-fP[1], m[2]-fP[2] };
1157
1158 {
1159 fvec mAVi[6] = { fC[0]-mV[0], fC[1]-mV[1], fC[2]-mV[2],
1160 fC[3]-mV[3], fC[4]-mV[4], fC[5]-mV[5] };
1161
1162// fvec init = InvertSym3( mAVi, mAVi ) ;
1163 InvertCholetsky3(mAVi);
1164
1165 fvec dChi2 = ( +(mAVi[0]*z[0] + mAVi[1]*z[1] + mAVi[3]*z[2])*z[0]
1166 +(mAVi[1]*z[0] + mAVi[2]*z[1] + mAVi[4]*z[2])*z[1]
1167 +(mAVi[3]*z[0] + mAVi[4]*z[1] + mAVi[5]*z[2])*z[2] );
1168
1169 // Take Abs(dChi2) here. Negative value of 'det' or 'dChi2' shows that the particle
1170 // was not used in the production vertex fit
1171
1172// fChi2 += (!init) & fabs( dChi2 );
1173 fChi2 += fabs( dChi2 );
1174 fNDF += 2;
1175 }
1176
1177 fP[0] = m[0];
1178 fP[1] = m[1];
1179 fP[2] = m[2];
1180 fP[3]+= mB[0][0]*z[0] + mB[0][1]*z[1] + mB[0][2]*z[2];
1181 fP[4]+= mB[1][0]*z[0] + mB[1][1]*z[1] + mB[1][2]*z[2];
1182 fP[5]+= mB[2][0]*z[0] + mB[2][1]*z[1] + mB[2][2]*z[2];
1183 fP[6]+= mB[3][0]*z[0] + mB[3][1]*z[1] + mB[3][2]*z[2];
1184 fP[7]+= mB[4][0]*z[0] + mB[4][1]*z[1] + mB[4][2]*z[2];
1185
1186 fvec d0, d1, d2;
1187
1188 fC[0] = mV[0];
1189 fC[1] = mV[1];
1190 fC[2] = mV[2];
1191 fC[3] = mV[3];
1192 fC[4] = mV[4];
1193 fC[5] = mV[5];
1194
1195 d0= mB[0][0]*mV[0] + mB[0][1]*mV[1] + mB[0][2]*mV[3] - fC[ 6];
1196 d1= mB[0][0]*mV[1] + mB[0][1]*mV[2] + mB[0][2]*mV[4] - fC[ 7];
1197 d2= mB[0][0]*mV[3] + mB[0][1]*mV[4] + mB[0][2]*mV[5] - fC[ 8];
1198
1199 fC[ 6]+= d0;
1200 fC[ 7]+= d1;
1201 fC[ 8]+= d2;
1202 fC[ 9]+= d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
1203
1204 d0= mB[1][0]*mV[0] + mB[1][1]*mV[1] + mB[1][2]*mV[3] - fC[10];
1205 d1= mB[1][0]*mV[1] + mB[1][1]*mV[2] + mB[1][2]*mV[4] - fC[11];
1206 d2= mB[1][0]*mV[3] + mB[1][1]*mV[4] + mB[1][2]*mV[5] - fC[12];
1207
1208 fC[10]+= d0;
1209 fC[11]+= d1;
1210 fC[12]+= d2;
1211 fC[13]+= d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
1212 fC[14]+= d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
1213
1214 d0= mB[2][0]*mV[0] + mB[2][1]*mV[1] + mB[2][2]*mV[3] - fC[15];
1215 d1= mB[2][0]*mV[1] + mB[2][1]*mV[2] + mB[2][2]*mV[4] - fC[16];
1216 d2= mB[2][0]*mV[3] + mB[2][1]*mV[4] + mB[2][2]*mV[5] - fC[17];
1217
1218 fC[15]+= d0;
1219 fC[16]+= d1;
1220 fC[17]+= d2;
1221 fC[18]+= d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
1222 fC[19]+= d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
1223 fC[20]+= d0*mB[2][0] + d1*mB[2][1] + d2*mB[2][2];
1224
1225 d0= mB[3][0]*mV[0] + mB[3][1]*mV[1] + mB[3][2]*mV[3] - fC[21];
1226 d1= mB[3][0]*mV[1] + mB[3][1]*mV[2] + mB[3][2]*mV[4] - fC[22];
1227 d2= mB[3][0]*mV[3] + mB[3][1]*mV[4] + mB[3][2]*mV[5] - fC[23];
1228
1229 fC[21]+= d0;
1230 fC[22]+= d1;
1231 fC[23]+= d2;
1232 fC[24]+= d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
1233 fC[25]+= d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
1234 fC[26]+= d0*mB[2][0] + d1*mB[2][1] + d2*mB[2][2];
1235 fC[27]+= d0*mB[3][0] + d1*mB[3][1] + d2*mB[3][2];
1236
1237 d0= mB[4][0]*mV[0] + mB[4][1]*mV[1] + mB[4][2]*mV[3] - fC[28];
1238 d1= mB[4][0]*mV[1] + mB[4][1]*mV[2] + mB[4][2]*mV[4] - fC[29];
1239 d2= mB[4][0]*mV[3] + mB[4][1]*mV[4] + mB[4][2]*mV[5] - fC[30];
1240
1241 fC[28]+= d0;
1242 fC[29]+= d1;
1243 fC[30]+= d2;
1244 fC[31]+= d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
1245 fC[32]+= d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
1246 fC[33]+= d0*mB[2][0] + d1*mB[2][1] + d2*mB[2][2];
1247 fC[34]+= d0*mB[3][0] + d1*mB[3][1] + d2*mB[3][2];
1248 fC[35]+= d0*mB[4][0] + d1*mB[4][1] + d2*mB[4][2];
1249
1250 if( noS ){
1251 fP[7] = 0;
1252 fC[28] = fC[29] = fC[30] = fC[31] = fC[32] = fC[33] = fC[34] = fC[35] = 0;
1253 } else {
1254 TransportToDS( fP[7] );
1255 Convert(0);
1256 }
1257
1258 fSFromDecay = 0;
1259}
1260
1261void KFParticleBaseSIMD::SetMassConstraint( fvec *mP, fvec *mC, fvec mJ[7][7], fvec mass, fvec mask )
1262{
1263 //* Set nonlinear mass constraint (Mass) on the state vector mP with a covariance matrix mC.
1264
1265 const fvec energy2 = mP[6]*mP[6], p2 = mP[3]*mP[3]+mP[4]*mP[4]+mP[5]*mP[5], mass2 = mass*mass;
1266
1267 const fvec a = energy2 - p2 + 2.*mass2;
1268 const fvec b = -2.*(energy2 + p2);
1269 const fvec c = energy2 - p2 - mass2;
1270
1271 fvec lambda = if3( fvec( fabs(b) > fvec(1.e-10) ), -c/b , Zero);
1272
1273 fvec d = 4.*energy2*p2 - mass2*(energy2-p2-2.*mass2);
1274 fvec qMask = fvec( fvec(d >= Zero) & fvec(fabs(a) > fvec(1.e-10)) );
1275 lambda = if3( qMask, ((energy2 + p2 - sqrt(d))/a), lambda );
1276
1277 lambda = if3( fvec(mP[6]<Zero), fvec(-1000000.f), lambda );
1278
1279 Int_t iIter=0;
1280 for(iIter=0; iIter<100; iIter++)
1281 {
1282 fvec lambda2 = lambda*lambda;
1283 fvec lambda4 = lambda2*lambda2;
1284
1285// fvec lambda0 = lambda;
1286
1287 fvec f = -mass2 * lambda4 + a*lambda2 + b*lambda + c;
1288 fvec df = -4.*mass2 * lambda2*lambda + 2.*a*lambda + b;
1289 lambda -= if3(fvec( fabs(df) > fvec(1.e-10) ) , f/df, Zero);
1290// if(TMath::Abs(lambda0 - lambda) < 1.e-8) break;
1291 }
1292
1293 const fvec lpi = 1./(1. + lambda);
1294 const fvec lmi = 1./(1. - lambda);
1295 const fvec lp2i = lpi*lpi;
1296 const fvec lm2i = lmi*lmi;
1297
1298 fvec lambda2 = lambda*lambda;
1299
1300 fvec dfl = -4.*mass2 * lambda2*lambda + 2.*a*lambda + b;
1301 fvec dfx[7] = {0};//,0,0,0};
1302 dfx[0] = -2.*(1. + lambda)*(1. + lambda)*mP[3];
1303 dfx[1] = -2.*(1. + lambda)*(1. + lambda)*mP[4];
1304 dfx[2] = -2.*(1. + lambda)*(1. + lambda)*mP[5];
1305 dfx[3] = 2.*(1. - lambda)*(1. - lambda)*mP[6];
1306 fvec dlx[4] = {1,1,1,1};
1307
1308 for(int i=0; i<4; i++)
1309 dlx[i] = if3( fvec( fabs(dfl) > fvec(1.e-10) ), -dfx[i] / dfl, dlx[i] );
1310
1311 fvec dxx[4] = {mP[3]*lm2i, mP[4]*lm2i, mP[5]*lm2i, -mP[6]*lp2i};
1312
1313 for(Int_t i=0; i<7; i++)
1314 for(Int_t j=0; j<7; j++)
1315 mJ[i][j]=0;
1316 mJ[0][0] = 1.;
1317 mJ[1][1] = 1.;
1318 mJ[2][2] = 1.;
1319
1320 for(Int_t i=3; i<7; i++)
1321 for(Int_t j=3; j<7; j++)
1322 mJ[i][j] = dlx[j-3]*dxx[i-3];
1323
1324 for(Int_t i=3; i<6; i++)
1325 mJ[i][i] += lmi;
1326 mJ[6][6] += lpi;
1327
1328 fvec mCJ[7][7];
1329
1330 for(Int_t i=0; i<7; i++) {
1331 for(Int_t j=0; j<7; j++) {
1332 mCJ[i][j] = 0;
1333 for(Int_t k=0; k<7; k++) {
1334 mCJ[i][j] += mC[IJ(i,k)]*mJ[j][k];
1335 }
1336 }
1337 }
1338
1339 for(Int_t i=0; i<7; ++i){
1340 for(Int_t j=0; j<=i; ++j){
1341 mC[IJ(i,j)] = if3(mask, Zero, mC[IJ(i,j)]);
1342 for(Int_t l=0; l<7; l++){
1343 mC[IJ(i,j)] += if3(mask, mJ[i][l]*mCJ[l][j], Zero);
1344 }
1345 }
1346 }
1347
1348 mP[3] *= if3(mask, lmi, One);
1349 mP[4] *= if3(mask, lmi, One);
1350 mP[5] *= if3(mask, lmi, One);
1351 mP[6] *= if3(mask, lpi, One);
1352}
1353
1355{
1356 //* Set nonlinear mass constraint (mass)
1357
1358 fvec mJ[7][7];
1359 SetMassConstraint( fP, fC, mJ, mass, fvec( One > Zero) );
1360 fMassHypo = mass;
1361 SumDaughterMass = mass;
1362}
1363
1365{
1366 //* Set hard( SigmaMass=0 ) or soft (SigmaMass>0) mass constraint
1367
1368 fMassHypo = Mass;
1369 SumDaughterMass = Mass;
1370
1371 fvec m2 = Mass*Mass; // measurement, weighted by Mass
1372 fvec s2 = m2*SigmaMass*SigmaMass; // sigma^2
1373
1374 fvec p2 = fP[3]*fP[3] + fP[4]*fP[4] + fP[5]*fP[5];
1375 fvec e0 = sqrt(m2+p2);
1376
1377 fvec mH[8];
1378 mH[0] = mH[1] = mH[2] = 0.;
1379 mH[3] = -2*fP[3];
1380 mH[4] = -2*fP[4];
1381 mH[5] = -2*fP[5];
1382 mH[6] = 2*fP[6];//e0;
1383 mH[7] = 0;
1384
1385 fvec zeta = e0*e0 - e0*fP[6];
1386 zeta = m2 - (fP[6]*fP[6]-p2);
1387
1388 fvec mCHt[8], s2_est=0;
1389 for( Int_t i=0; i<8; ++i ){
1390 mCHt[i] = 0.0;
1391 for (Int_t j=0;j<8;++j) mCHt[i] += Cij(i,j)*mH[j];
1392 s2_est += mH[i]*mCHt[i];
1393 }
1394
1395//TODO Next line should be uncommented
1396// if( s2_est<1.e-20 ) return; // calculated mass error is already 0,
1397 // the particle can not be constrained on mass
1398
1399 fvec w2 = 1./( s2 + s2_est );
1400 fChi2 += zeta*zeta*w2;
1401 fNDF += 1;
1402 for( Int_t i=0, ii=0; i<8; ++i ){
1403 fvec ki = mCHt[i]*w2;
1404 fP[i]+= ki*zeta;
1405 for(Int_t j=0;j<=i;++j) fC[ii++] -= ki*mCHt[j];
1406 }
1407}
1408
1409
1411{
1412 //* Set no decay length for resonances
1413
1415
1416 fvec h[8];
1417 h[0] = h[1] = h[2] = h[3] = h[4] = h[5] = h[6] = 0;
1418 h[7] = 1;
1419
1420 fvec zeta = 0 - fP[7];
1421 for(Int_t i=0;i<8;++i) zeta -= h[i]*(fP[i]-fP[i]);
1422
1423 fvec s = fC[35];
1424// if( s>1.e-20 ) //TODO Should be uncommented!
1425 {
1426 s = 1./s;
1427 fChi2 += zeta*zeta*s;
1428 fNDF += 1;
1429 for( Int_t i=0, ii=0; i<7; ++i ){
1430 fvec ki = fC[28+i]*s;
1431 fP[i]+= ki*zeta;
1432 for(Int_t j=0;j<=i;++j) fC[ii++] -= ki*fC[28+j];
1433 }
1434 }
1435 fP[7] = 0;
1436 fC[28] = fC[29] = fC[30] = fC[31] = fC[32] = fC[33] = fC[34] = fC[35] = 0;
1437}
1438
1439void KFParticleBaseSIMD::Construct( const KFParticleBaseSIMD* vDaughters[], Int_t nDaughters,
1440 const KFParticleBaseSIMD *Parent, Float_t Mass, Bool_t IsConstrained ,
1441 Bool_t isAtVtxGuess )
1442{
1443 //* Full reconstruction in one go
1444
1445 //* Parameters: array of pointers to the input daughter particles (vDaughters), a number of daughters (nDaughters),
1446 //* pointer to the production vertex or particle (Parent), if the pointer is set - a constraint on the production
1447 //* is set on the reconstructed particle, if 0 - will not be set, a mass hypothesis (Mass), if Mass>0 - a mass constraint
1448 //* will be set, a flag (IsConstrained), which indicates whether daughters are constrained on the production region,
1449 //* a flag (isAtVtxGuess), which indicates whether daughters are already transported to the vertex guess or not,
1450 //* if it is set to 1, a vertex guess should be provided, by default is set to 0.
1451
1452 //* examples:
1453 //* 1) construction of a secondary particle with no vertex guess:
1454 //* const KFParticleSIMD* vDaughtersPointer[2] = {&vDaughters[0], &vDaughters[1]};
1455 //* mother.Construct(vDaughtersPointer, 2, 0);
1456 //* 2) construction of the primary particle with a vertex guess
1457 //* fvec errGuess[3] = {100*sqrt(PrimVtx.CovarianceMatrix()[0]),
1458 //* 100*sqrt(PrimVtx.CovarianceMatrix()[2]),
1459 //* 100*sqrt(PrimVtx.CovarianceMatrix()[5])};
1460 //* mother.SetVtxGuess(PrimVtx.X(), PrimVtx.Y(), PrimVtx.Z());
1461 //* mother.SetVtxErrGuess(errGuess[0], errGuess[1], errGuess[2]);
1462 //* const KFParticleSIMD* vDaughtersPointer[2] = {&vDaughters[0], &vDaughters[1]};
1463 //* mother.Construct(vDaughtersPointer, 2, 0, -1, 0, 1);
1464
1465 Int_t maxIter = 1;
1466 bool wasLinearized = fIsLinearized;
1467 if( (!fIsLinearized || IsConstrained) && (!isAtVtxGuess) ){
1468 //fVtxGuess[0] = fVtxGuess[1] = fVtxGuess[2] = 0; //!!!!
1469 fvec ds=0., ds1=0.;
1470 fvec P[8], C[36];
1471 vDaughters[0]->GetDStoParticle(*vDaughters[1], ds, ds1);
1472 vDaughters[0]->Transport(ds,P,C);
1473 fVtxGuess[0] = P[0];
1474 fVtxGuess[1] = P[1];
1475 fVtxGuess[2] = P[2];
1476/* fVtxGuess[0] = GetX();
1477 fVtxGuess[1] = GetY();
1478 fVtxGuess[2] = GetZ();*/
1479
1480 if(!fIsVtxErrGuess)
1481 {
1482 fVtxErrGuess[0] = if3(fvec(C[0]>Zero), 10*sqrt(C[0]), One);
1483 fVtxErrGuess[1] = if3(fvec(C[2]>Zero), 10*sqrt(C[2]), One);
1484 fVtxErrGuess[2] = if3(fvec(C[5]>Zero), 10*sqrt(C[5]), One);
1485 }
1486
1487 fIsLinearized = 1;
1488 maxIter = 3;
1489 }
1490 else {
1491 if(!fIsVtxErrGuess)
1492 {
1493 fVtxErrGuess[0] = 1;
1494 fVtxErrGuess[1] = 1;
1495 fVtxErrGuess[2] = 1;
1496 }
1497 }
1498
1499 fvec constraintC[6];
1500
1501 if( IsConstrained ){
1502 for(Int_t i=0;i<6;++i) constraintC[i]=fC[i];
1503 }
1504 else {
1505 for(Int_t i=0;i<6;++i) constraintC[i]=0.;
1506// constraintC[0] = constraintC[2] = constraintC[5] = 100.;
1507
1508 constraintC[0] = fVtxErrGuess[0]*fVtxErrGuess[0];
1509 constraintC[2] = fVtxErrGuess[1]*fVtxErrGuess[1];
1510 constraintC[5] = fVtxErrGuess[2]*fVtxErrGuess[2];
1511 }
1512
1513
1514 for( Int_t iter=0; iter<maxIter; iter++ ){
1515
1517 SetNDaughters(nDaughters);
1518
1520 fSFromDecay = 0;
1521 fP[0] = fVtxGuess[0];
1522 fP[1] = fVtxGuess[1];
1523 fP[2] = fVtxGuess[2];
1524 fP[3] = 0;
1525 fP[4] = 0;
1526 fP[5] = 0;
1527 fP[6] = 0;
1528 fP[7] = 0;
1529 SumDaughterMass = 0;
1530
1531 for(Int_t i=0;i<6; ++i) fC[i]=constraintC[i];
1532 for(Int_t i=6;i<36;++i) fC[i]=0.;
1533 fC[35] = 1.;
1534
1535 fNDF = IsConstrained ?0 :-3;
1536 fChi2 = 0.;
1537 fQ = 0;
1538
1539 for( Int_t itr =0; itr<nDaughters; itr++ ){
1540 AddDaughter( *vDaughters[itr], isAtVtxGuess );
1541 }
1542 if( iter<maxIter-1){
1543 for( Int_t i=0; i<3; i++ ) fVtxGuess[i] = fP[i];
1544 }
1545 }
1546 fIsLinearized = wasLinearized;
1547
1548 if( Mass>=0 ) SetMassConstraint( Mass );
1549 if( Parent ) SetProductionVertex( *Parent );
1550}
1551
1552void KFParticleBaseSIMD::Convert( bool ToProduction )
1553{
1554 //* Tricky function - convert the particle error along its trajectory to
1555 //* the value which corresponds to its production/decay vertex
1556 //* It is done by combination of the error of decay length with the position errors
1557
1558 fvec fld[3];
1559 {
1560 GetFieldValue( fP, fld );
1561 const fvec kCLight = fQ*0.000299792458;
1562 fld[0]*=kCLight; fld[1]*=kCLight; fld[2]*=kCLight;
1563 }
1564
1565 fvec h[6];
1566
1567 h[0] = fP[3];
1568 h[1] = fP[4];
1569 h[2] = fP[5];
1570 if( ToProduction ){ h[0]=-h[0]; h[1]=-h[1]; h[2]=-h[2]; }
1571 h[3] = h[1]*fld[2]-h[2]*fld[1];
1572 h[4] = h[2]*fld[0]-h[0]*fld[2];
1573 h[5] = h[0]*fld[1]-h[1]*fld[0];
1574
1575 fvec c;
1576
1577 c = fC[28]+h[0]*fC[35];
1578 fC[ 0]+= h[0]*(c+fC[28]);
1579 fC[28] = c;
1580
1581 fC[ 1]+= h[1]*fC[28] + h[0]*fC[29];
1582 c = fC[29]+h[1]*fC[35];
1583 fC[ 2]+= h[1]*(c+fC[29]);
1584 fC[29] = c;
1585
1586 fC[ 3]+= h[2]*fC[28] + h[0]*fC[30];
1587 fC[ 4]+= h[2]*fC[29] + h[1]*fC[30];
1588 c = fC[30]+h[2]*fC[35];
1589 fC[ 5]+= h[2]*(c+fC[30]);
1590 fC[30] = c;
1591
1592 fC[ 6]+= h[3]*fC[28] + h[0]*fC[31];
1593 fC[ 7]+= h[3]*fC[29] + h[1]*fC[31];
1594 fC[ 8]+= h[3]*fC[30] + h[2]*fC[31];
1595 c = fC[31]+h[3]*fC[35];
1596 fC[ 9]+= h[3]*(c+fC[31]);
1597 fC[31] = c;
1598
1599 fC[10]+= h[4]*fC[28] + h[0]*fC[32];
1600 fC[11]+= h[4]*fC[29] + h[1]*fC[32];
1601 fC[12]+= h[4]*fC[30] + h[2]*fC[32];
1602 fC[13]+= h[4]*fC[31] + h[3]*fC[32];
1603 c = fC[32]+h[4]*fC[35];
1604 fC[14]+= h[4]*(c+fC[32]);
1605 fC[32] = c;
1606
1607 fC[15]+= h[5]*fC[28] + h[0]*fC[33];
1608 fC[16]+= h[5]*fC[29] + h[1]*fC[33];
1609 fC[17]+= h[5]*fC[30] + h[2]*fC[33];
1610 fC[18]+= h[5]*fC[31] + h[3]*fC[33];
1611 fC[19]+= h[5]*fC[32] + h[4]*fC[33];
1612 c = fC[33]+h[5]*fC[35];
1613 fC[20]+= h[5]*(c+fC[33]);
1614 fC[33] = c;
1615
1616 fC[21]+= h[0]*fC[34];
1617 fC[22]+= h[1]*fC[34];
1618 fC[23]+= h[2]*fC[34];
1619 fC[24]+= h[3]*fC[34];
1620 fC[25]+= h[4]*fC[34];
1621 fC[26]+= h[5]*fC[34];
1622}
1623
1624
1626{
1627 //* Transport the particle to its decay vertex
1628
1632}
1633
1635{
1636 //* Transport the particle to its production vertex
1637
1639 if( !fAtProductionVertex ) Convert( 1 );
1641}
1642
1643
1645{
1646 //* Transport the particle on dS parameter (SignedPath/Momentum)
1647
1648 Transport( dS, fP, fC );
1649 fSFromDecay+= dS;
1650}
1651
1652
1653void KFParticleBaseSIMD::GetDistanceToVertexLine( const KFParticleBaseSIMD &Vertex, fvec &l, fvec &dl, fvec *isParticleFromVertex ) const
1654{
1655 //* Get dS to a certain space point without field
1656
1657 fvec c[6] = {Vertex.fC[0]+fC[0], Vertex.fC[1]+fC[1], Vertex.fC[2]+fC[2],
1658 Vertex.fC[3]+fC[3], Vertex.fC[4]+fC[4], Vertex.fC[5]+fC[5]};
1659
1660 fvec dx = (Vertex.fP[0]-fP[0]);
1661 fvec dy = (Vertex.fP[1]-fP[1]);
1662 fvec dz = (Vertex.fP[2]-fP[2]);
1663
1664
1665 l = sqrt( dx*dx + dy*dy + dz*dz );
1666 dl = c[0]*dx*dx + c[2]*dy*dy + c[5]*dz*dz + 2*(c[1]*dx*dy + c[3]*dx*dz + c[4]*dy*dz);
1667 fvec ok = fvec(fvec(0)<dl);
1668 dl = fvec(ok & (sqrt( dl )/l)) + fvec(!ok & fvec(1.e8));
1669
1670 if(isParticleFromVertex)
1671 {
1672 *isParticleFromVertex = fvec(ok & fvec( l<fvec(3*dl) ));
1673 fvec cos = dx*fP[3] + dy*fP[4] + dz*fP[5];
1674// fvec dCos = dy*dy*fC[14] + dz*dz*fC[20] + dx*dx*fC[9] + 2*dz*fC[15]*fP[3] + c[0]* fP[3]*fP[3] +
1675// 2*dz*fC[16]* fP[4] + 2 *c[1] *fP[3] *fP[4] + c[2] *fP[4]*fP[4] + 2 *dz *fC[17]* fP[5] +
1676// 2*c[3] *fP[3]* fP[5] + 2 *c[4] *fP[4] *fP[5] + c[5]*fP[5] *fP[5] +
1677// 2*dy *(dz *fC[19] + fC[10] *fP[3] + fC[11]* fP[4] + fC[12]* fP[5]) +
1678// 2*dx *(dy *fC[13] + dz *fC[18] + fC[6]* fP[3] + fC[7]* fP[4] + fC[8]* fP[5]);
1679// ok = fvec(fvec(0)<dCos);
1680// dCos = fvec(ok & ( dCos ));
1681// dCos = sqrt(dCos);
1682 *isParticleFromVertex = fvec( (*isParticleFromVertex) | fvec(!(*isParticleFromVertex) & fvec(cos<Zero) ) );
1683 }
1684}
1685
1686
1688 const
1689{
1690
1691 //* Get dS to a certain space point for Bz field
1692 const fvec kCLight = 0.000299792458;
1693 fvec bq = B*fQ*kCLight;
1694 fvec pt2 = fP[3]*fP[3] + fP[4]*fP[4];
1695
1696 fvec dx = xyz[0] - fP[0];
1697 fvec dy = xyz[1] - fP[1];
1698 fvec a = dx*fP[3]+dy*fP[4];
1699 fvec dS = Zero;
1700
1701 const fvec LocalSmall = 1.e-8;
1702 fvec mask = fvec( fabs(bq)<LocalSmall );
1703 dS = if3(mask, a/pt2, (atan2( bq*a, pt2 + bq*(dy*fP[3] -dx*fP[4]) )/bq));
1704
1705 #if 0
1706 {
1707
1708 fvec px = fP[3];
1709 fvec py = fP[4];
1710 fvec pz = fP[5];
1711 fvec ss[2], g[2][5];
1712
1713 ss[0] = dS;
1714 ss[1] = -dS;
1715 for( Int_t i=0; i<2; i++){
1716 fvec bs = bq*ss[i];
1717 fvec c = TMath::Cos(bs), s = TMath::Sin(bs);
1718 fvec cB,sB;
1719 if( TMath::Abs(bq)>1.e-8){
1720 cB= (1-c)/bq;
1721 sB= s/bq;
1722 }else{
1723 const fvec kOvSqr6 = 1./TMath::Sqrt(6.);
1724 sB = (1.-bs*kOvSqr6)*(1.+bs*kOvSqr6)*ss[i];
1725 cB = .5*sB*bs;
1726 }
1727 g[i][0] = fP[0] + sB*px + cB*py;
1728 g[i][1] = fP[1] - cB*px + sB*py;
1729 g[i][2] = fP[2] + ss[i]*pz;
1730 g[i][3] = + c*px + s*py;
1731 g[i][4] = - s*px + c*py;
1732 }
1733
1734 Int_t i=0;
1735
1736 fvec dMin = 1.e10;
1737 for( Int_t j=0; j<2; j++){
1738 fvec xx = g[j][0]-xyz[0];
1739 fvec yy = g[j][1]-xyz[1];
1740 fvec zz = g[j][2]-xyz[2];
1741 fvec d = xx*xx + yy*yy + zz*zz;
1742 if( d<dMin ){
1743 dMin = d;
1744 i = j;
1745 }
1746 }
1747
1748 dS = ss[i];
1749
1750 fvec x= g[i][0], y= g[i][1], z= g[i][2], ppx= g[i][3], ppy= g[i][4];
1751 fvec ddx = x-xyz[0];
1752 fvec ddy = y-xyz[1];
1753 fvec ddz = z-xyz[2];
1754 fvec c = ddx*ppx + ddy*ppy + ddz*pz ;
1755 fvec pp2 = ppx*ppx + ppy*ppy + pz*pz;
1756 if( TMath::Abs(pp2)>1.e-8 ){
1757 dS+=c/pp2;
1758 }
1759 }
1760 #endif
1761 return dS;
1762}
1763
1765 const
1766{
1767
1768 //* Get dS to a certain space point for Bz field
1769
1770 const fvec kCLight = 0.000299792458;
1771 fvec bq = B*fQ*kCLight;
1772 fvec pt2 = fP[3]*fP[3] + fP[5]*fP[5];
1773
1774 fvec dx = xyz[0] - fP[0];
1775 fvec dy = -xyz[2] + fP[2];
1776 fvec a = dx*fP[3]-dy*fP[5];
1777 fvec dS = Zero;
1778
1779 const fvec LocalSmall = 1.e-8;
1780 fvec mask = fvec( fabs(bq)<LocalSmall );
1781 dS = if3(mask, a/pt2, (atan2( bq*a, pt2 + bq*(dy*fP[3] +dx*fP[5]) )/bq));
1782 return dS;
1783}
1784
1786{
1787 //* Get dS to another particle for Bz field
1788
1789 fvec px = fP[3];
1790 fvec py = fP[4];
1791 fvec pz = fP[5];
1792
1793 fvec px1 = p.fP[3];
1794 fvec py1 = p.fP[4];
1795 fvec pz1 = p.fP[5];
1796
1797 const fvec kCLight = 0.000299792458;
1798
1799 fvec bq = B*fQ*kCLight;
1800 fvec bq1 = B*p.fQ*kCLight;
1801 fvec s=Zero, ds=Zero, s1=Zero, ds1=Zero;
1802
1803 const fvec LocalSmall = 1.e-8;
1804
1805 fvec bq_big = fvec(LocalSmall < fabs(bq));
1806 fvec bq1_big = fvec(LocalSmall < fabs(bq1));
1807 {
1808 const fvec two = 2.;
1809 const fvec four = 2.;
1810 fvec dx = (p.fP[0] - fP[0]);
1811 fvec dy = (p.fP[1] - fP[1]);
1812 fvec d2 = (dx*dx+dy*dy);
1813
1814 fvec p2 = (px *px + py *py);
1815 fvec p21 = (px1*px1 + py1*py1);
1816
1817 fvec a = (px*py1 - py*px1);
1818 fvec b = (px*px1 + py*py1);
1819
1820 fvec ldx = bq*bq1*dx - bq1*py + bq*py1 ;
1821 fvec ldy = bq*bq1*dy + bq1*px - bq*px1 ;
1822 fvec l2 = ldx*ldx + ldy*ldy;
1823
1824 fvec cS = bq1*p2 + bq*bq1*(dy* px - dx* py) - bq*b;
1825 fvec cS1= bq*p21 - bq*bq1*(dy*px1 - dx*py1) - bq1*b;
1826
1827 fvec ca = bq*bq*bq1*d2 + two*( cS + bq*bq*(py1*dx-px1*dy)) ;
1828 fvec ca1 = bq*bq1*bq1*d2 + two*( cS1 - bq1*bq1*(py*dx-px*dy)) ;
1829
1830 fvec sa = 4*l2*p2 - ca*ca;
1831 fvec sa1 = 4*l2*p21 - ca1*ca1;
1832
1833 const fvec saNeg = fvec(sa<Zero);
1834 sa = (!saNeg) & sa;
1835 const fvec sa1Neg = fvec(sa1<Zero);
1836 sa1= (!sa1Neg) & sa1;
1837
1838
1839 s = if3( bq_big, (atan2( bq*( bq1*(dx*px +dy*py) + a ) , cS )/bq), (( (dx*px + dy*py) + (py*px1-px*py1)/bq1)/p2) );
1840 ds = if3( bq_big, (atan2(sqrt(sa),ca)/bq) , (s*s - (d2-two*(px1*dy-py1*dx)/bq1)/p2) );
1841 fvec dsNeg = fvec(ds<Zero);
1842 ds = ( fvec((!dsNeg) & (!bq_big)) & sqrt(ds) ) + (bq_big & ds);
1843
1844 s1 = if3( bq1_big, (atan2( -bq1*( bq*(dx*px1+dy*py1) + a), cS1 )/bq1), ((-(dx*px1 + dy*py1) + (py*px1-px*py1)/bq)/p21) );
1845 ds1 = if3( bq1_big, (atan2(sqrt(sa1),ca1)/bq1) , (s1*s1 - (d2+two*(px*dy-py*dx)/bq)/p21) );
1846 fvec ds1Neg = fvec(ds1<Zero);
1847 ds1 = ( fvec((!ds1Neg) & (!bq1_big)) & sqrt(ds1) ) + (bq1_big & ds1);
1848 }
1849
1850 fvec ss[2], ss1[2], g[2][5],g1[2][5];
1851
1852 ss[0] = s + ds;
1853 ss[1] = s - ds;
1854 ss1[0] = s1 + ds1;
1855 ss1[1] = s1 - ds1;
1856// std::cout << " ss " << ss[0][0] << " " << ss[1][0] << " " << ss1[0][0] << " " << ss1[1][0] << std::endl;
1857 for( Int_t i=0; i<2; i++){
1858 fvec bs = bq*ss[i];
1859 fvec c = cos(bs), sss = sin(bs);
1860
1861 fvec cB,sB;
1862 const fvec kOvSqr6 = 1./sqrt(6.);
1863 sB = if3( bq_big, (sss/bq), ((One-bs*kOvSqr6)*(One+bs*kOvSqr6)*ss[i]));
1864 cB = if3( bq_big, ((One-c)/bq), (.5*sB*bs));
1865
1866 g[i][0] = fP[0] + sB*px + cB*py;
1867 g[i][1] = fP[1] - cB*px + sB*py;
1868 g[i][2] = fP[2] + ss[i]*pz;
1869 g[i][3] = + c*px + sss*py;
1870 g[i][4] = - sss*px + c*py;
1871
1872 bs = bq1*ss1[i];
1873 c = cos(bs); sss = sin(bs);
1874
1875 sB = if3( bq1_big, (sss/bq1), ((One-bs*kOvSqr6)*(One+bs*kOvSqr6)*ss1[i]));
1876 cB = if3( bq1_big, ((One-c)/bq1), (.5*sB*bs));
1877
1878 g1[i][0] = p.fP[0] + sB*px1 + cB*py1;
1879 g1[i][1] = p.fP[1] - cB*px1 + sB*py1;
1880 g1[i][2] = p.fP[2] + ss[i]*pz1;
1881 g1[i][3] = + c*px1 + sss*py1;
1882 g1[i][4] = - sss*px1 + c*py1;
1883 }
1884
1885 DS = ss[0];
1886 DS1 = ss1[0];
1887
1888 fvec dMin = 1.e10;
1889 for( Int_t j=0; j<2; j++){
1890 for( Int_t j1=0; j1<2; j1++){
1891 fvec xx = g[j][0]-g1[j1][0];
1892 fvec yy = g[j][1]-g1[j1][1];
1893 fvec zz = g[j][2]-g1[j1][2];
1894 fvec d = xx*xx + yy*yy + zz*zz;
1895
1896 fvec mask = fvec(d<dMin);
1897 DS = (mask & ss[j]) + ((!mask) & DS);
1898 DS1 = (mask & ss1[j]) + ((!mask) & DS1);
1899 dMin = (mask & d) + ((!mask) & dMin);
1900 }
1901 }
1902 #if 0
1903 {
1904 fvec x= g[i][0], y= g[i][1], z= g[i][2], ppx= g[i][3], ppy= g[i][4];
1905 fvec x1=g1[i1][0], y1= g1[i1][1], z1= g1[i1][2], ppx1= g1[i1][3], ppy1= g1[i1][4];
1906 fvec dx = x1-x;
1907 fvec dy = y1-y;
1908 fvec dz = z1-z;
1909 fvec a = ppx*ppx1 + ppy*ppy1 + pz*pz1;
1910 fvec b = dx*ppx1 + dy*ppy1 + dz*pz1;
1911 fvec c = dx*ppx + dy*ppy + dz*pz ;
1912 fvec pp2 = ppx*ppx + ppy*ppy + pz*pz;
1913 fvec pp21= ppx1*ppx1 + ppy1*ppy1 + pz1*pz1;
1914 fvec det = pp2*pp21 - a*a;
1915 if( fabs(det)>1.e-8 ){
1916 DS+=(a*b-pp21*c)/det;
1917 DS1+=(a*c-pp2*b)/det;
1918 }
1919 }
1920 #endif
1921}
1922
1924{
1925 //* Get dS to another particle for Bz field
1926
1927 fvec px = fP[3];
1928 fvec py = -fP[5];
1929 fvec pz = fP[4];
1930
1931 fvec px1 = p.fP[3];
1932 fvec py1 = -p.fP[5];
1933 fvec pz1 = p.fP[4];
1934
1935 const fvec kCLight = 0.000299792458;
1936
1937 fvec bq = B*fQ*kCLight;
1938 fvec bq1 = B*p.fQ*kCLight;
1939 fvec s=Zero, ds=Zero, s1=Zero, ds1=Zero;
1940
1941 const fvec LocalSmall = 1.e-8;
1942
1943 fvec bq_big = fvec(LocalSmall < fabs(bq));
1944 fvec bq1_big = fvec(LocalSmall < fabs(bq1));
1945 {
1946 const fvec two = 2.;
1947 const fvec four = 2.;
1948 fvec dx = (p.fP[0] - fP[0]);
1949 fvec dy = (-p.fP[2] + fP[2]);
1950 fvec d2 = (dx*dx+dy*dy);
1951
1952 fvec p2 = (px *px + py *py);
1953 fvec p21 = (px1*px1 + py1*py1);
1954
1955 fvec a = (px*py1 - py*px1);
1956 fvec b = (px*px1 + py*py1);
1957
1958 fvec ldx = bq*bq1*dx - bq1*py + bq*py1 ;
1959 fvec ldy = bq*bq1*dy + bq1*px - bq*px1 ;
1960 fvec l2 = ldx*ldx + ldy*ldy;
1961
1962 fvec cS = bq1*p2 + bq*bq1*(dy* px - dx* py) - bq*b;
1963 fvec cS1= bq*p21 - bq*bq1*(dy*px1 - dx*py1) - bq1*b;
1964
1965 fvec ca = bq*bq*bq1*d2 + two*( cS + bq*bq*(py1*dx-px1*dy)) ;
1966 fvec ca1 = bq*bq1*bq1*d2 + two*( cS1 - bq1*bq1*(py*dx-px*dy)) ;
1967
1968 fvec sa = 4*l2*p2 - ca*ca;
1969 fvec sa1 = 4*l2*p21 - ca1*ca1;
1970
1971 const fvec saNeg = fvec(sa<Zero);
1972 sa = (!saNeg) & sa;
1973 const fvec sa1Neg = fvec(sa1<Zero);
1974 sa1= (!sa1Neg) & sa1;
1975
1976 s = if3( bq_big, (atan2( bq*( bq1*(dx*px +dy*py) + a ) , cS )/bq), (( (dx*px + dy*py) + (py*px1-px*py1)/bq1)/p2) );
1977 ds = if3( bq_big, (atan2(sqrt(sa),ca)/bq) , (s*s - (d2-two*(px1*dy-py1*dx)/bq1)/p2) );
1978 fvec dsNeg = fvec(ds<Zero);
1979 ds = ( fvec((!dsNeg) & (!bq_big)) & sqrt(ds) ) + (bq_big & ds);
1980
1981 s1 = if3( bq1_big, (atan2( -bq1*( bq*(dx*px1+dy*py1) + a), cS1 )/bq1), ((-(dx*px1 + dy*py1) + (py*px1-px*py1)/bq)/p21) );
1982 ds1 = if3( bq1_big, (atan2(sqrt(sa1),ca1)/bq1) , (s1*s1 - (d2+two*(px*dy-py*dx)/bq)/p21) );
1983 fvec ds1Neg = fvec(ds1<Zero);
1984 ds1 = ( fvec((!ds1Neg) & (!bq1_big)) & sqrt(ds1) ) + (bq1_big & ds1);
1985 }
1986
1987 fvec ss[2], ss1[2], g[2][5],g1[2][5];
1988
1989 ss[0] = s + ds;
1990 ss[1] = s - ds;
1991 ss1[0] = s1 + ds1;
1992 ss1[1] = s1 - ds1;
1993// std::cout << " ss " << ss[0][0] << " " << ss[1][0] << " " << ss1[0][0] << " " << ss1[1][0] << std::endl;
1994 for( Int_t i=0; i<2; i++){
1995 fvec bs = bq*ss[i];
1996 fvec c = cos(bs), sss = sin(bs);
1997
1998 fvec cB,sB;
1999 const fvec kOvSqr6 = 1./sqrt(6.);
2000 sB = if3( bq_big, (sss/bq), ((One-bs*kOvSqr6)*(One+bs*kOvSqr6)*ss[i]));
2001 cB = if3( bq_big, ((One-c)/bq), (.5*sB*bs));
2002
2003 g[i][0] = fP[0] + sB*px + cB*py;
2004 g[i][1] = -fP[2] - cB*px + sB*py;
2005 g[i][2] = fP[1] + ss[i]*pz;
2006 g[i][3] = + c*px + sss*py;
2007 g[i][4] = - sss*px + c*py;
2008
2009 bs = bq1*ss1[i];
2010 c = cos(bs); sss = sin(bs);
2011
2012 sB = if3( bq1_big, (sss/bq1), ((One-bs*kOvSqr6)*(One+bs*kOvSqr6)*ss1[i]));
2013 cB = if3( bq1_big, ((One-c)/bq1), (.5*sB*bs));
2014
2015 g1[i][0] = p.fP[0] + sB*px1 + cB*py1;
2016 g1[i][1] = -p.fP[2] - cB*px1 + sB*py1;
2017 g1[i][2] = p.fP[1] + ss[i]*pz1;
2018 g1[i][3] = + c*px1 + sss*py1;
2019 g1[i][4] = - sss*px1 + c*py1;
2020 }
2021
2022 DS = ss[0];
2023 DS1 = ss1[0];
2024
2025 fvec dMin = 1.e10;
2026 for( Int_t j=0; j<2; j++){
2027 for( Int_t j1=0; j1<2; j1++){
2028 fvec xx = g[j][0]-g1[j1][0];
2029 fvec yy = g[j][1]-g1[j1][1];
2030 fvec zz = g[j][2]-g1[j1][2];
2031 fvec d = xx*xx + yy*yy + zz*zz;
2032
2033 fvec mask = fvec(d<dMin);
2034 DS = (mask & ss[j]) + ((!mask) & DS);
2035 DS1 = (mask & ss1[j]) + ((!mask) & DS1);
2036 dMin = (mask & d) + ((!mask) & dMin);
2037 }
2038 }
2039}
2040
2042{
2043 //* Transport the particle on dS, output to P[],C[], for CBM field
2044
2045 const fvec *r = fP;
2046 fvec p2 = r[3]*r[3] + r[4]*r[4] + r[5]*r[5];
2047 fvec dS = if3(fabs(p2) > fvec(1.E-4), (( r[3]*(xyz[0]-r[0]) + r[4]*(xyz[1]-r[1]) + r[5]*(xyz[2]-r[2]) )/p2), Zero);
2048
2049// fvec dS = 0;
2050//
2051// if( fQ[0]==0 ){
2052// dS = GetDStoPointLine( xyz );
2053// return dS;
2054// }
2055//
2056// fvec fld[3];
2057// GetFieldValue( fP, fld );
2058// dS = GetDStoPointBy( fld[1],xyz );
2059
2060 dS = if3(fabs(dS)>1.E3, Zero, dS);
2061
2062 return dS;
2063}
2064
2066{
2067 fvec p12 = fP[3]*fP[3] + fP[4]*fP[4] + fP[5]*fP[5];
2068 fvec p22 = p.fP[3]*p.fP[3] + p.fP[4]*p.fP[4] + p.fP[5]*p.fP[5];
2069 fvec p1p2 = fP[3]*p.fP[3] + fP[4]*p.fP[4] + fP[5]*p.fP[5];
2070
2071 fvec drp1 = fP[3]*(p.fP[0]-fP[0]) + fP[4]*(p.fP[1]-fP[1]) + fP[5]*(p.fP[2]-fP[2]);
2072 fvec drp2 = p.fP[3]*(p.fP[0]-fP[0]) + p.fP[4]*(p.fP[1]-fP[1]) + p.fP[5]*(p.fP[2]-fP[2]);
2073
2074 fvec detp = p1p2*p1p2 - p12*p22;
2075
2076 dS = (drp2*p1p2 - drp1*p22) /detp; //TODO put defence against detp = 0
2077 dS1 = (drp2*p12 - drp1*p1p2)/detp;
2078}
2079
2081{
2082 //* Transport the particle on dS, output to P[],C[], for CBM field
2083
2084 fvec mTy[2];
2085 fvec mY[2];
2086 fvec mZ[2];
2087
2088 mY[0] = GetY();
2089 mZ[0] = GetZ();
2090 mTy[0] = if3( fabs(GetPz()) > small, GetPy()/GetPz(), 1.);
2091
2092 mY[1] = p.GetY();
2093 mZ[1] = p.GetZ();
2094 mTy[1] = if3( fabs(p.GetPz()) > small, p.GetPy()/p.GetPz(), 1.);
2095
2096 fvec r0[3];
2097 fvec dty = mTy[0]-mTy[1];
2098
2099 r0[2] = if3( fabs(dty) > small, (mTy[0]*mZ[0]-mTy[1]*mZ[1] + mY[1] -mY[0])/dty, 0. );
2100 r0[1] = mY[0] + mTy[0]*(r0[2]-mZ[0]);
2101 r0[0] = GetX() + GetPx()/GetPz()*(r0[2]-mZ[0]);
2102
2103 dS = GetDStoPointCBM(r0);
2104 dS1 = p.GetDStoPoint(r0);
2105
2106// if( fQ[0]==0 ){
2107// GetDStoParticleLine( p, dS, dS1 );
2108// return;
2109// }
2110//
2111// fvec fld[3];
2112// GetFieldValue( fP, fld );
2113//
2114// GetDStoParticleBy(fld[1],p,dS,dS1);
2115}
2116
2118 fvec P[], fvec C[] ) const
2119{
2120 //* Transport the particle on dS, output to P[],C[], for CBM field
2121
2122 if( fQ[0]==0 ){
2123 TransportLine( dS, P, C );
2124 return;
2125 }
2126
2127 const fvec kCLight = 0.000299792458;
2128
2129 fvec c = fQ*kCLight;
2130
2131 // construct coefficients
2132
2133 fvec
2134 px = fP[3],
2135 py = fP[4],
2136 pz = fP[5];
2137
2138 fvec sx=0, sy=0, sz=0, syy=0, syz=0, syyy=0, ssx=0, ssy=0, ssz=0, ssyy=0, ssyz=0, ssyyy=0;
2139
2140 { // get field integrals
2141
2142 fvec fld[3][3];
2143 fvec p0[3], p1[3], p2[3];
2144
2145 // line track approximation
2146
2147 p0[0] = fP[0];
2148 p0[1] = fP[1];
2149 p0[2] = fP[2];
2150
2151 p2[0] = fP[0] + px*dS;
2152 p2[1] = fP[1] + py*dS;
2153 p2[2] = fP[2] + pz*dS;
2154
2155 p1[0] = 0.5*(p0[0]+p2[0]);
2156 p1[1] = 0.5*(p0[1]+p2[1]);
2157 p1[2] = 0.5*(p0[2]+p2[2]);
2158
2159 // first order track approximation
2160 {
2161 GetFieldValue( p0, fld[0] );
2162 GetFieldValue( p1, fld[1] );
2163 GetFieldValue( p2, fld[2] );
2164
2165 fvec ssy1 = ( 7*fld[0][1] + 6*fld[1][1]-fld[2][1] )*c*dS*dS/96.;
2166 fvec ssy2 = ( fld[0][1] + 2*fld[1][1] )*c*dS*dS/6.;
2167
2168 p1[0] -= ssy1*pz;
2169 p1[2] += ssy1*px;
2170 p2[0] -= ssy2*pz;
2171 p2[2] += ssy2*px;
2172 }
2173
2174 GetFieldValue( p0, fld[0] );
2175 GetFieldValue( p1, fld[1] );
2176 GetFieldValue( p2, fld[2] );
2177
2178 for(int iF1=0; iF1<3; iF1++)
2179 for(int iF2=0; iF2<3; iF2++)
2180 fld[iF1][iF2] = if3( (fabs(fld[iF1][iF2]) > fvec(100.)), Zero, fld[iF1][iF2] );
2181
2182 sx = c*( fld[0][0] + 4*fld[1][0] + fld[2][0] )*dS/6.;
2183 sy = c*( fld[0][1] + 4*fld[1][1] + fld[2][1] )*dS/6.;
2184 sz = c*( fld[0][2] + 4*fld[1][2] + fld[2][2] )*dS/6.;
2185
2186 ssx = c*( fld[0][0] + 2*fld[1][0])*dS*dS/6.;
2187 ssy = c*( fld[0][1] + 2*fld[1][1])*dS*dS/6.;
2188 ssz = c*( fld[0][2] + 2*fld[1][2])*dS*dS/6.;
2189
2190 fvec c2[3][3] = { { 5, -4, -1},{ 44, 80, -4},{ 11, 44, 5} }; // /=360.
2191 fvec cc2[3][3] = { { 38, 8, -4},{ 148, 208, -20},{ 3, 36, 3} }; // /=2520.
2192 for(Int_t n=0; n<3; n++)
2193 for(Int_t m=0; m<3; m++)
2194 {
2195 syz += c2[n][m]*fld[n][1]*fld[m][2];
2196 ssyz += cc2[n][m]*fld[n][1]*fld[m][2];
2197 }
2198
2199 syz *= c*c*dS*dS/360.;
2200 ssyz *= c*c*dS*dS*dS/2520.;
2201
2202 syy = c*( fld[0][1] + 4*fld[1][1] + fld[2][1] )*dS;
2203 syyy = syy*syy*syy / 1296;
2204 syy = syy*syy/72;
2205
2206 ssyy = ( fld[0][1]*( 38*fld[0][1] + 156*fld[1][1] - fld[2][1] )+
2207 fld[1][1]*( 208*fld[1][1] +16*fld[2][1] )+
2208 fld[2][1]*( 3*fld[2][1] )
2209 )*dS*dS*dS*c*c/2520.;
2210 ssyyy =
2211 (
2212 fld[0][1]*( fld[0][1]*( 85*fld[0][1] + 526*fld[1][1] - 7*fld[2][1] )+
2213 fld[1][1]*( 1376*fld[1][1] +84*fld[2][1] )+
2214 fld[2][1]*( 19*fld[2][1] ) )+
2215 fld[1][1]*( fld[1][1]*( 1376*fld[1][1] +256*fld[2][1] )+
2216 fld[2][1]*( 62*fld[2][1] ) )+
2217 fld[2][1]*fld[2][1] *( 3*fld[2][1] )
2218 )*dS*dS*dS*dS*c*c*c/90720.;
2219
2220 }
2221
2222 fvec mJ[11];
2223
2224 mJ[0]=dS-ssyy; mJ[1]=ssx; mJ[2]=ssyyy-ssy;
2225 mJ[3]=-ssz; mJ[4]=dS; mJ[5]=ssx+ssyz;
2226
2227 mJ[6]=1-syy; mJ[7]=sx; mJ[8]=syyy-sy;
2228 mJ[9]=-sz; mJ[10]=sx+syz;
2229
2230
2231
2232 P[0] = fP[0] + mJ[0]*px + mJ[1]*py + mJ[2]*pz;
2233 P[1] = fP[1] + mJ[3]*px + mJ[4]*py + mJ[5]*pz;
2234 P[2] = fP[2] - mJ[2]*px - mJ[1]*py + mJ[0]*pz;
2235 P[3] = mJ[6]*px + mJ[7]*py + mJ[8]*pz;
2236 P[4] = mJ[9]*px + py + mJ[10]*pz;
2237 P[5] = -mJ[8]*px - mJ[7]*py + mJ[6]*pz;
2238 P[6] = fP[6];
2239 P[7] = fP[7];
2240
2241 if(C!=fC)
2242 {
2243 for(int iC=0; iC<36; iC++)
2244 C[iC] = fC[iC];
2245 }
2246
2247 multQSQt1( mJ, C);
2248
2249// fvec mJ[8][8];
2250// for( Int_t i=0; i<8; i++ ) for( Int_t j=0; j<8; j++) mJ[i][j]=0;
2251//
2252// 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;
2253// 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;
2254// 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;
2255//
2256// 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;
2257// 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;
2258// 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;
2259// mJ[6][6] = mJ[7][7] = 1;
2260//
2261// P[0] = fP[0] + mJ[0][3]*px + mJ[0][4]*py + mJ[0][5]*pz;
2262// P[1] = fP[1] + mJ[1][3]*px + mJ[1][4]*py + mJ[1][5]*pz;
2263// P[2] = fP[2] + mJ[2][3]*px + mJ[2][4]*py + mJ[2][5]*pz;
2264// P[3] = mJ[3][3]*px + mJ[3][4]*py + mJ[3][5]*pz;
2265// P[4] = mJ[4][3]*px + mJ[4][4]*py + mJ[4][5]*pz;
2266// P[5] = mJ[5][3]*px + mJ[5][4]*py + mJ[5][5]*pz;
2267// P[6] = fP[6];
2268// P[7] = fP[7];
2269
2270// MultQSQt( mJ, fC, C);
2271
2272}
2273
2275{
2276 const fvec A00 = S[0 ]+S[6 ]*J[0]+S[10]*J[1]+S[15]*J[2];
2277 const fvec A10 = S[1 ]+S[7 ]*J[0]+S[11]*J[1]+S[16]*J[2];
2278 const fvec A20 = S[3 ]+S[8 ]*J[0]+S[12]*J[1]+S[17]*J[2];
2279 const fvec A30 = S[6 ]+S[9 ]*J[0]+S[13]*J[1]+S[18]*J[2];
2280 const fvec A40 = S[10]+S[13]*J[0]+S[14]*J[1]+S[19]*J[2];
2281 const fvec A50 = S[15]+S[18]*J[0]+S[19]*J[1]+S[20]*J[2];
2282 const fvec A60 = S[21]+S[24]*J[0]+S[25]*J[1]+S[26]*J[2];
2283 const fvec A70 = S[28]+S[31]*J[0]+S[32]*J[1]+S[33]*J[2];
2284
2285 S[0 ] = A00+J[0]*A30+J[1]*A40+J[2]*A50;
2286 S[1 ] = A10+J[3]*A30+J[4]*A40+J[5]*A50;
2287 S[3 ] = A20-J[2]*A30-J[1]*A40+J[0]*A50;
2288 S[6 ] = J[6]*A30+J[7]*A40+J[8]*A50;
2289 S[10] = J[9]*A30+ A40+J[10]*A50;
2290 S[15] = -J[8]*A30-J[7]*A40+J[6]*A50;
2291 S[21] = A60;
2292 S[28] = A70;
2293
2294// const fvec A01 = S[1 ]+S[6 ]*J[3]+S[10]*J[4]+S[15]*J[5];
2295 const fvec A11 = S[2 ]+S[7 ]*J[3]+S[11]*J[4]+S[16]*J[5];
2296 const fvec A21 = S[4 ]+S[8 ]*J[3]+S[12]*J[4]+S[17]*J[5];
2297 const fvec A31 = S[7 ]+S[9 ]*J[3]+S[13]*J[4]+S[18]*J[5];
2298 const fvec A41 = S[11]+S[13]*J[3]+S[14]*J[4]+S[19]*J[5];
2299 const fvec A51 = S[16]+S[18]*J[3]+S[19]*J[4]+S[20]*J[5];
2300 const fvec A61 = S[22]+S[24]*J[3]+S[25]*J[4]+S[26]*J[5];
2301 const fvec A71 = S[29]+S[31]*J[3]+S[32]*J[4]+S[33]*J[5];
2302
2303 S[2 ] = A11+J[3]*A31+J[4]*A41+J[5]*A51;
2304 S[4 ] = A21-J[2]*A31-J[1]*A41+J[0]*A51;
2305 S[7 ] = J[6]*A31+J[7]*A41+J[8]*A51;
2306 S[11] = J[9]*A31+ A41+J[10]*A51;
2307 S[16] = -J[8]*A31-J[7]*A41+J[6]*A51;
2308 S[22] = A61;
2309 S[29] = A71;
2310
2311// const fvec A02 = S[3 ]+S[6 ]*J[2][3]+S[10]*J[2][4]+S[15]*J[2][5];
2312// const fvec A12 = S[4 ]+S[7 ]*J[2][3]+S[11]*J[2][4]+S[16]*J[2][5];
2313 const fvec A22 = S[5 ]-S[8 ]*J[2]-S[12]*J[1]+S[17]*J[0];
2314 const fvec A32 = S[8 ]-S[9 ]*J[2]-S[13]*J[1]+S[18]*J[0];
2315 const fvec A42 = S[12]-S[13]*J[2]-S[14]*J[1]+S[19]*J[0];
2316 const fvec A52 = S[17]-S[18]*J[2]-S[19]*J[1]+S[20]*J[0];
2317 const fvec A62 = S[23]-S[24]*J[2]-S[25]*J[1]+S[26]*J[0];
2318 const fvec A72 = S[30]-S[31]*J[2]-S[32]*J[1]+S[33]*J[0];
2319
2320 S[5 ] = A22-J[2]*A32-J[1]*A42+J[0]*A52;
2321 S[8 ] = J[6]*A32+J[7]*A42+J[8]*A52;
2322 S[12] = J[9]*A32+ A42+J[10]*A52;
2323 S[17] = -J[8]*A32-J[7]*A42+J[6]*A52;
2324 S[23] = A62;
2325 S[30] = A72;
2326
2327// const fvec A03 = S[6] *J[6]+S[10]*J[7]+S[15]*J[8];
2328// const fvec A13 = S[7] *J[6]+S[11]*J[7]+S[16]*J[8];
2329// const fvec A23 = S[8] *J[6]+S[12]*J[7]+S[17]*J[8];
2330 const fvec A33 = S[9] *J[6]+S[13]*J[7]+S[18]*J[8];
2331 const fvec A43 = S[13]*J[6]+S[14]*J[7]+S[19]*J[8];
2332 const fvec A53 = S[18]*J[6]+S[19]*J[7]+S[20]*J[8];
2333 const fvec A63 = S[24]*J[6]+S[25]*J[7]+S[26]*J[8];
2334 const fvec A73 = S[31]*J[6]+S[32]*J[7]+S[33]*J[8];
2335
2336// const fvec A04 = S[6] *J[9]+S[10]*J[4][4]+S[15]*J[10];
2337// const fvec A14 = S[7] *J[9]+S[11]*J[4][4]+S[16]*J[10];
2338// const fvec A24 = S[8] *J[9]+S[12]*J[4][4]+S[17]*J[10];
2339 const fvec A34 = S[9] *J[9]+S[13]+S[18]*J[10];
2340 const fvec A44 = S[13]*J[9]+S[14]+S[19]*J[10];
2341 const fvec A54 = S[18]*J[9]+S[19]+S[20]*J[10];
2342 const fvec A64 = S[24]*J[9]+S[25]+S[26]*J[10];
2343 const fvec A74 = S[31]*J[9]+S[32]+S[33]*J[10];
2344
2345// const fvec A05 = S[6] *J[5][3]+S[10]*J[5][4]+S[15]*J[6];
2346// const fvec A15 = S[7] *J[5][3]+S[11]*J[5][4]+S[16]*J[6];
2347// const fvec A25 = S[8] *J[5][3]+S[12]*J[5][4]+S[17]*J[6];
2348 const fvec A35 = -S[9] *J[8]-S[13]*J[7]+S[18]*J[6];
2349 const fvec A45 = -S[13]*J[8]-S[14]*J[7]+S[19]*J[6];
2350 const fvec A55 = -S[18]*J[8]-S[19]*J[7]+S[20]*J[6];
2351 const fvec A65 = -S[24]*J[8]-S[25]*J[7]+S[26]*J[6];
2352 const fvec A75 = -S[31]*J[8]-S[32]*J[7]+S[33]*J[6];
2353
2354
2355 S[9 ] = J[6]*A33+J[7]*A43+J[8]*A53;
2356 S[13] = J[9]*A33+ A43+J[10]*A53;
2357 S[18] = -J[8]*A33-J[7]*A43+J[6]*A53;
2358 S[24] = A63;
2359 S[31] = A73;
2360
2361 S[14] = J[9]*A34+ A44+J[10]*A54;
2362 S[19] = -J[8]*A34-J[7]*A44+J[6]*A54;
2363 S[25] = A64;
2364 S[32] = A74;
2365
2366 S[20] = -J[8]*A35-J[7]*A45+J[6]*A55;
2367 S[26] = A65;
2368 S[33] = A75;
2369
2370//S[27] = S27;
2371//S[34] = S34;
2372//S[35] = S35;
2373}
2374
2376 fvec p[], fvec e[] ) const
2377{
2378 //* Transport the particle on dS, output to P[],C[], for Bz field
2379
2380 const fvec kCLight = 0.000299792458;
2381 b = b*fQ*kCLight;
2382 fvec bs= b*t;
2383 fvec s = sin(bs), c = cos(bs);
2384 fvec sB, cB;
2385
2386 const fvec kOvSqr6 = 1./sqrt(6.);
2387 const fvec LocalSmall = 1.e-10;
2388
2389 sB = if3( fvec(LocalSmall < fabs(bs)), (s/b), ((One-bs*kOvSqr6)*(One+bs*kOvSqr6)*t) );
2390 cB = if3( fvec(LocalSmall < fabs(bs)), ((One-c)/b), (.5*sB*bs) );
2391
2392 fvec px = fP[3];
2393 fvec py = fP[4];
2394 fvec pz = fP[5];
2395
2396 p[0] = fP[0] + sB*px + cB*py;
2397 p[1] = fP[1] - cB*px + sB*py;
2398 p[2] = fP[2] + t*pz;
2399 p[3] = c*px + s*py;
2400 p[4] = -s*px + c*py;
2401 p[5] = fP[5];
2402 p[6] = fP[6];
2403 p[7] = fP[7];
2404
2405 /*
2406 fvec mJ[8][8] = { {1,0,0, sB, cB, 0, 0, 0 },
2407 {0,1,0, -cB, sB, 0, 0, 0 },
2408 {0,0,1, 0, 0, t, 0, 0 },
2409 {0,0,0, c, s, 0, 0, 0 },
2410 {0,0,0, -s, c, 0, 0, 0 },
2411 {0,0,0, 0, 0, 1, 0, 0 },
2412 {0,0,0, 0, 0, 0, 1, 0 },
2413 {0,0,0, 0, 0, 0, 0, 1 } };
2414 fvec mA[8][8];
2415 for( Int_t k=0,i=0; i<8; i++)
2416 for( Int_t j=0; j<=i; j++, k++ ) mA[i][j] = mA[j][i] = fC[k];
2417
2418 fvec mJC[8][8];
2419 for( Int_t i=0; i<8; i++ )
2420 for( Int_t j=0; j<8; j++ ){
2421 mJC[i][j]=0;
2422 for( Int_t k=0; k<8; k++ ) mJC[i][j]+=mJ[i][k]*mA[k][j];
2423 }
2424
2425 for( Int_t k=0,i=0; i<8; i++)
2426 for( Int_t j=0; j<=i; j++, k++ ){
2427 e[k] = 0;
2428 for( Int_t l=0; l<8; l++ ) e[k]+=mJC[i][l]*mJ[j][l];
2429 }
2430
2431 return;
2432 */
2433
2434 fvec
2435 c6=fC[6], c7=fC[7], c8=fC[8], c17=fC[17], c18=fC[18],
2436 c24 = fC[24], c31 = fC[31];
2437
2438 fvec
2439 cBC13 = cB*fC[13],
2440 mJC13 = c7 - cB*fC[9] + sB*fC[13],
2441 mJC14 = fC[11] - cBC13 + sB*fC[14],
2442 mJC23 = c8 + t*c18,
2443 mJC24 = fC[12] + t*fC[19],
2444 mJC33 = c*fC[9] + s*fC[13],
2445 mJC34 = c*fC[13] + s*fC[14],
2446 mJC43 = -s*fC[9] + c*fC[13],
2447 mJC44 = -s*fC[13] + c*fC[14];
2448
2449
2450 e[0]= fC[0] + 2*(sB*c6 + cB*fC[10]) + (sB*fC[9] + 2*cBC13)*sB + cB*cB*fC[14];
2451 e[1]= fC[1] - cB*c6 + sB*fC[10] + mJC13*sB + mJC14*cB;
2452 e[2]= fC[2] - cB*c7 + sB*fC[11] - mJC13*cB + mJC14*sB;
2453 e[3]= fC[3] + t*fC[15] + mJC23*sB + mJC24*cB;
2454 e[4]= fC[4] + t*fC[16] - mJC23*cB + mJC24*sB;
2455
2456 e[15]= fC[15] + c18*sB + fC[19]*cB;
2457 e[16]= fC[16] - c18*cB + fC[19]*sB;
2458 e[17]= c17 + fC[20]*t;
2459 e[18]= c18*c + fC[19]*s;
2460 e[19]= -c18*s + fC[19]*c;
2461
2462 e[5]= fC[5] + (c17 + e[17] )*t;
2463
2464 e[6]= c*c6 + s*fC[10] + mJC33*sB + mJC34*cB;
2465 e[7]= c*c7 + s*fC[11] - mJC33*cB + mJC34*sB;
2466 e[8]= c*c8 + s*fC[12] + e[18]*t;
2467 e[9]= mJC33*c + mJC34*s;
2468 e[10]= -s*c6 + c*fC[10] + mJC43*sB + mJC44*cB;
2469
2470
2471 e[11]= -s*c7 + c*fC[11] - mJC43*cB + mJC44*sB;
2472 e[12]= -s*c8 + c*fC[12] + e[19]*t;
2473 e[13]= mJC43*c + mJC44*s;
2474 e[14]= -mJC43*s + mJC44*c;
2475 e[20]= fC[20];
2476 e[21]= fC[21] + fC[25]*cB + c24*sB;
2477 e[22]= fC[22] - c24*cB + fC[25]*sB;
2478 e[23]= fC[23] + fC[26]*t;
2479 e[24]= c*c24 + s*fC[25];
2480 e[25]= c*fC[25] - c24*s;
2481 e[26]= fC[26];
2482 e[27]= fC[27];
2483 e[28]= fC[28] + fC[32]*cB + c31*sB;
2484 e[29]= fC[29] - c31*cB + fC[32]*sB;
2485 e[30]= fC[30] + fC[33]*t;
2486 e[31]= c*c31 + s*fC[32];
2487 e[32]= c*fC[32] - s*c31;
2488 e[33]= fC[33];
2489 e[34]= fC[34];
2490 e[35]= fC[35];
2491}
2492
2493
2495{
2496 //* Calculate distance from vertex [cm]
2497
2498 return GetDistanceFromVertex( Vtx.fP );
2499}
2500
2502{
2503 //* Calculate distance from vertex [cm]
2504
2505 fvec mP[8], mC[36];
2506 Transport( GetDStoPoint(vtx), mP, mC );
2507 fvec d[3]={ vtx[0]-mP[0], vtx[1]-mP[1], vtx[2]-mP[2]};
2508 return sqrt( d[0]*d[0]+d[1]*d[1]+d[2]*d[2] );
2509}
2510
2512 const
2513{
2514 //* Calculate distance to other particle [cm]
2515
2516 fvec dS, dS1;
2517 GetDStoParticle( p, dS, dS1 );
2518 fvec mP[8], mC[36], mP1[8], mC1[36];
2519 Transport( dS, mP, mC );
2520 p.Transport( dS1, mP1, mC1 );
2521 fvec dx = mP[0]-mP1[0];
2522 fvec dy = mP[1]-mP1[1];
2523 fvec dz = mP[2]-mP1[2];
2524 return sqrt(dx*dx+dy*dy+dz*dz);
2525}
2526
2528{
2529 //* Calculate sqrt(Chi2/ndf) deviation from vertex
2530
2531 return GetDeviationFromVertex( Vtx.fP, Vtx.fC );
2532}
2533
2534
2536{
2537 //* Calculate sqrt(Chi2/ndf) deviation from vertex
2538 //* v = [xyz], Cv=[Cxx,Cxy,Cyy,Cxz,Cyz,Czz]-covariance matrix
2539
2540 fvec mP[8];
2541 fvec mC[36];
2542
2543 Transport( GetDStoPoint(v), mP, mC );
2544
2545 fvec d[3]={ v[0]-mP[0], v[1]-mP[1], v[2]-mP[2]};
2546
2547 fvec sigmaS = .1+10.*sqrt( (d[0]*d[0]+d[1]*d[1]+d[2]*d[2])/
2548 (mP[3]*mP[3]+mP[4]*mP[4]+mP[5]*mP[5]) );
2549
2550
2551 fvec h[3] = { mP[3]*sigmaS, mP[4]*sigmaS, mP[5]*sigmaS };
2552
2553 fvec mSi[6] =
2554 { mC[0] +h[0]*h[0],
2555 mC[1] +h[1]*h[0], mC[2] +h[1]*h[1],
2556 mC[3] +h[2]*h[0], mC[4] +h[2]*h[1], mC[5] +h[2]*h[2] };
2557
2558 if( Cv ){
2559 mSi[0]+=Cv[0];
2560 mSi[1]+=Cv[1];
2561 mSi[2]+=Cv[2];
2562 mSi[3]+=Cv[3];
2563 mSi[4]+=Cv[4];
2564 mSi[5]+=Cv[5];
2565 }
2566
2567 fvec mS[6];
2568
2569 mS[0] = mSi[2]*mSi[5] - mSi[4]*mSi[4];
2570 mS[1] = mSi[3]*mSi[4] - mSi[1]*mSi[5];
2571 mS[2] = mSi[0]*mSi[5] - mSi[3]*mSi[3];
2572 mS[3] = mSi[1]*mSi[4] - mSi[2]*mSi[3];
2573 mS[4] = mSi[1]*mSi[3] - mSi[0]*mSi[4];
2574 mS[5] = mSi[0]*mSi[2] - mSi[1]*mSi[1];
2575
2576 fvec s = ( mSi[0]*mS[0] + mSi[1]*mS[1] + mSi[3]*mS[3] );
2577 s = if3( fvec( fabs(s) > 1.E-8 ) , One/s , Zero);
2578
2579 return sqrt( fabs(s*( ( mS[0]*d[0] + mS[1]*d[1] + mS[3]*d[2])*d[0]
2580 +(mS[1]*d[0] + mS[2]*d[1] + mS[4]*d[2])*d[1]
2581 +(mS[3]*d[0] + mS[4]*d[1] + mS[5]*d[2])*d[2] ))/2);
2582}
2583
2584
2586 const
2587{
2588 //* Calculate sqrt(Chi2/ndf) deviation from other particle
2589
2590 fvec dS, dS1;
2591 GetDStoParticle( p, dS, dS1 );
2592 fvec mP1[8], mC1[36];
2593 p.Transport( dS1, mP1, mC1 );
2594
2595 fvec d[3]={ fP[0]-mP1[0], fP[1]-mP1[1], fP[2]-mP1[2]};
2596
2597 fvec sigmaS = .1+10.*sqrt( (d[0]*d[0]+d[1]*d[1]+d[2]*d[2])/
2598 (mP1[3]*mP1[3]+mP1[4]*mP1[4]+mP1[5]*mP1[5]) );
2599
2600 fvec h[3] = { mP1[3]*sigmaS, mP1[4]*sigmaS, mP1[5]*sigmaS };
2601
2602 mC1[0] +=h[0]*h[0];
2603 mC1[1] +=h[1]*h[0];
2604 mC1[2] +=h[1]*h[1];
2605 mC1[3] +=h[2]*h[0];
2606 mC1[4] +=h[2]*h[1];
2607 mC1[5] +=h[2]*h[2];
2608
2609 return GetDeviationFromVertex( mP1, mC1 )*sqrt(2./1.);
2610}
2611
2612
2613
2615{
2616 //* Subtract the particle from the vertex
2617
2618 fvec m[8];
2619 fvec mCm[36];
2620
2621 if( Vtx.fIsLinearized ){
2622 GetMeasurement( Vtx.fVtxGuess, m, mCm );
2623 } else {
2624 GetMeasurement( Vtx.fP, m, mCm );
2625 }
2626
2627 fvec mV[6];
2628
2629 mV[ 0] = mCm[ 0];
2630 mV[ 1] = mCm[ 1];
2631 mV[ 2] = mCm[ 2];
2632 mV[ 3] = mCm[ 3];
2633 mV[ 4] = mCm[ 4];
2634 mV[ 5] = mCm[ 5];
2635
2636 //*
2637
2638 fvec mS[6];
2639 {
2640 fvec mSi[6] = { mV[0]-Vtx.fC[0],
2641 mV[1]-Vtx.fC[1], mV[2]-Vtx.fC[2],
2642 mV[3]-Vtx.fC[3], mV[4]-Vtx.fC[4], mV[5]-Vtx.fC[5] };
2643
2644 mS[0] = mSi[2]*mSi[5] - mSi[4]*mSi[4];
2645 mS[1] = mSi[3]*mSi[4] - mSi[1]*mSi[5];
2646 mS[2] = mSi[0]*mSi[5] - mSi[3]*mSi[3];
2647 mS[3] = mSi[1]*mSi[4] - mSi[2]*mSi[3];
2648 mS[4] = mSi[1]*mSi[3] - mSi[0]*mSi[4];
2649 mS[5] = mSi[0]*mSi[2] - mSi[1]*mSi[1];
2650
2651 fvec s = ( mSi[0]*mS[0] + mSi[1]*mS[1] + mSi[3]*mS[3] );
2652 s = if3(fvec( fabs(s) > small ) , One/s, Zero);
2653 mS[0]*=s;
2654 mS[1]*=s;
2655 mS[2]*=s;
2656 mS[3]*=s;
2657 mS[4]*=s;
2658 mS[5]*=s;
2659 }
2660
2661 //* Residual (measured - estimated)
2662
2663 fvec zeta[3] = { m[0]-Vtx.fP[0], m[1]-Vtx.fP[1], m[2]-Vtx.fP[2] };
2664
2665 //* mCHt = mCH' - D'
2666
2667 fvec mCHt0[3], mCHt1[3], mCHt2[3];
2668
2669 mCHt0[0]=Vtx.fC[ 0] ; mCHt1[0]=Vtx.fC[ 1] ; mCHt2[0]=Vtx.fC[ 3] ;
2670 mCHt0[1]=Vtx.fC[ 1] ; mCHt1[1]=Vtx.fC[ 2] ; mCHt2[1]=Vtx.fC[ 4] ;
2671 mCHt0[2]=Vtx.fC[ 3] ; mCHt1[2]=Vtx.fC[ 4] ; mCHt2[2]=Vtx.fC[ 5] ;
2672
2673 //* Kalman gain K = mCH'*S
2674
2675 fvec k0[3], k1[3], k2[3];
2676
2677 for(Int_t i=0;i<3;++i){
2678 k0[i] = mCHt0[i]*mS[0] + mCHt1[i]*mS[1] + mCHt2[i]*mS[3];
2679 k1[i] = mCHt0[i]*mS[1] + mCHt1[i]*mS[2] + mCHt2[i]*mS[4];
2680 k2[i] = mCHt0[i]*mS[3] + mCHt1[i]*mS[4] + mCHt2[i]*mS[5];
2681 }
2682
2683 //* New estimation of the vertex position r += K*zeta
2684
2685 fvec dChi2 = -(mS[0]*zeta[0] + mS[1]*zeta[1] + mS[3]*zeta[2])*zeta[0]
2686 - (mS[1]*zeta[0] + mS[2]*zeta[1] + mS[4]*zeta[2])*zeta[1]
2687 - (mS[3]*zeta[0] + mS[4]*zeta[1] + mS[5]*zeta[2])*zeta[2];
2688
2689// fvec mask = fvec( (Vtx.fChi2 - dChi2) < Zero ); //TODO Check and correct!
2690 fvec mask = fvec(One>Zero); //TODO Check and correct!
2691
2692 for(Int_t i=0;i<3;++i)
2693 Vtx.fP[i] -= (mask & (k0[i]*zeta[0] + k1[i]*zeta[1] + k2[i]*zeta[2]) );
2694 //* New covariance matrix C -= K*(mCH')'
2695
2696 for(Int_t i=0, k=0;i<3;++i){
2697 for(Int_t j=0;j<=i;++j,++k)
2698 Vtx.fC[k] += (mask & (k0[i]*mCHt0[j] + k1[i]*mCHt1[j] + k2[i]*mCHt2[j]));
2699 }
2700
2701 //* Calculate Chi^2
2702
2703 Vtx.fNDF -= 2; //TODO Check and correct!
2704 Vtx.fChi2 += (mask & dChi2);
2705}
2706
2708{
2709 //* Subtract the particle from the mother particle
2710
2711 fvec m[8];
2712 fvec mV[36];
2713
2714 if( Vtx.fIsLinearized ){
2715 GetMeasurement( Vtx.fVtxGuess, m, mV,0 );
2716 } else {
2717 GetMeasurement( Vtx.fP, m, mV );
2718 }
2719
2720 fvec mS[6]= { mV[0] - Vtx.fC[0],
2721 mV[1] - Vtx.fC[1], mV[2] - Vtx.fC[2],
2722 mV[3] - Vtx.fC[3], mV[4] - Vtx.fC[4], mV[5] - Vtx.fC[5] };
2723 InvertCholetsky3(mS);
2724
2725 //* Residual (measured - estimated)
2726
2727 fvec zeta[3] = { m[0]-Vtx.fP[0], m[1]-Vtx.fP[1], m[2]-Vtx.fP[2] };
2728
2729 //* CHt = CH' - D'
2730
2731 fvec mCHt0[7], mCHt1[7], mCHt2[7];
2732
2733 mCHt0[0]=mV[ 0] ; mCHt1[0]=mV[ 1] ; mCHt2[0]=mV[ 3] ;
2734 mCHt0[1]=mV[ 1] ; mCHt1[1]=mV[ 2] ; mCHt2[1]=mV[ 4] ;
2735 mCHt0[2]=mV[ 3] ; mCHt1[2]=mV[ 4] ; mCHt2[2]=mV[ 5] ;
2736 mCHt0[3]=Vtx.fC[ 6]-mV[ 6]; mCHt1[3]=Vtx.fC[ 7]-mV[ 7]; mCHt2[3]=Vtx.fC[ 8]-mV[ 8];
2737 mCHt0[4]=Vtx.fC[10]-mV[10]; mCHt1[4]=Vtx.fC[11]-mV[11]; mCHt2[4]=Vtx.fC[12]-mV[12];
2738 mCHt0[5]=Vtx.fC[15]-mV[15]; mCHt1[5]=Vtx.fC[16]-mV[16]; mCHt2[5]=Vtx.fC[17]-mV[17];
2739 mCHt0[6]=Vtx.fC[21]-mV[21]; mCHt1[6]=Vtx.fC[22]-mV[22]; mCHt2[6]=Vtx.fC[23]-mV[23];
2740
2741 //* Kalman gain K = mCH'*S
2742
2743 fvec k0[7], k1[7], k2[7];
2744
2745 for(Int_t i=0;i<7;++i){
2746 k0[i] = mCHt0[i]*mS[0] + mCHt1[i]*mS[1] + mCHt2[i]*mS[3];
2747 k1[i] = mCHt0[i]*mS[1] + mCHt1[i]*mS[2] + mCHt2[i]*mS[4];
2748 k2[i] = mCHt0[i]*mS[3] + mCHt1[i]*mS[4] + mCHt2[i]*mS[5];
2749 }
2750
2751 //* Add the daughter momentum to the particle momentum
2752
2753 Vtx.fP[ 3] -= m[ 3];
2754 Vtx.fP[ 4] -= m[ 4];
2755 Vtx.fP[ 5] -= m[ 5];
2756 Vtx.fP[ 6] -= m[ 6];
2757
2758 Vtx.fC[ 9] -= mV[ 9];
2759 Vtx.fC[13] -= mV[13];
2760 Vtx.fC[14] -= mV[14];
2761 Vtx.fC[18] -= mV[18];
2762 Vtx.fC[19] -= mV[19];
2763 Vtx.fC[20] -= mV[20];
2764 Vtx.fC[24] -= mV[24];
2765 Vtx.fC[25] -= mV[25];
2766 Vtx.fC[26] -= mV[26];
2767 Vtx.fC[27] -= mV[27];
2768
2769 //* New estimation of the vertex position r += K*zeta
2770
2771 for(Int_t i=0;i<3;++i)
2772 Vtx.fP[i] = m[i] - (k0[i]*zeta[0] + k1[i]*zeta[1] + k2[i]*zeta[2]);
2773 for(Int_t i=3;i<7;++i)
2774 Vtx.fP[i] = Vtx.fP[i] - (k0[i]*zeta[0] + k1[i]*zeta[1] + k2[i]*zeta[2]);
2775
2776 //* New covariance matrix C -= K*(mCH')'
2777
2778 fvec ffC[28] = { -mV[ 0],
2779 -mV[ 1], -mV[ 2],
2780 -mV[ 3], -mV[ 4], -mV[ 5],
2781 mV[ 6], mV[ 7], mV[ 8], Vtx.fC[ 9],
2782 mV[10], mV[11], mV[12], Vtx.fC[13], Vtx.fC[14],
2783 mV[15], mV[16], mV[17], Vtx.fC[18], Vtx.fC[19], Vtx.fC[20],
2784 mV[21], mV[22], mV[23], Vtx.fC[24], Vtx.fC[25], Vtx.fC[26], Vtx.fC[27] };
2785
2786 for(Int_t i=0, k=0;i<7;++i){
2787 for(Int_t j=0;j<=i;++j,++k){
2788 Vtx.fC[k] = ffC[k] + (k0[i]*mCHt0[j] + k1[i]*mCHt1[j] + k2[i]*mCHt2[j] );
2789 }
2790 }
2791
2792 //* Calculate Chi^2
2793 Vtx.fNDF -= 2;
2794 Vtx.fQ -= GetQ();
2795 Vtx.fSFromDecay = 0;
2796 Vtx.fChi2 -= ((mS[0]*zeta[0] + mS[1]*zeta[1] + mS[3]*zeta[2])*zeta[0]
2797 +(mS[1]*zeta[0] + mS[2]*zeta[1] + mS[4]*zeta[2])*zeta[1]
2798 +(mS[3]*zeta[0] + mS[4]*zeta[1] + mS[5]*zeta[2])*zeta[2]);
2799}
2800
2802 fvec P[], fvec C[] ) const
2803{
2804 //* Transport the particle as a straight line
2805
2806 P[0] = fP[0] + dS*fP[3];
2807 P[1] = fP[1] + dS*fP[4];
2808 P[2] = fP[2] + dS*fP[5];
2809 P[3] = fP[3];
2810 P[4] = fP[4];
2811 P[5] = fP[5];
2812 P[6] = fP[6];
2813 P[7] = fP[7];
2814
2815 fvec c6 = fC[ 6] + dS*fC[ 9];
2816 fvec c11 = fC[11] + dS*fC[14];
2817 fvec c17 = fC[17] + dS*fC[20];
2818 fvec sc13 = dS*fC[13];
2819 fvec sc18 = dS*fC[18];
2820 fvec sc19 = dS*fC[19];
2821
2822 C[ 0] = fC[ 0] + dS*( fC[ 6] + c6 );
2823 C[ 2] = fC[ 2] + dS*( fC[11] + c11 );
2824 C[ 5] = fC[ 5] + dS*( fC[17] + c17 );
2825
2826 C[ 7] = fC[ 7] + sc13;
2827 C[ 8] = fC[ 8] + sc18;
2828 C[ 9] = fC[ 9];
2829
2830 C[12] = fC[12] + sc19;
2831
2832 C[ 1] = fC[ 1] + dS*( fC[10] + C[ 7] );
2833 C[ 3] = fC[ 3] + dS*( fC[15] + C[ 8] );
2834 C[ 4] = fC[ 4] + dS*( fC[16] + C[12] );
2835 C[ 6] = c6;
2836
2837 C[10] = fC[10] + sc13;
2838 C[11] = c11;
2839
2840 C[13] = fC[13];
2841 C[14] = fC[14];
2842 C[15] = fC[15] + sc18;
2843 C[16] = fC[16] + sc19;
2844 C[17] = c17;
2845
2846 C[18] = fC[18];
2847 C[19] = fC[19];
2848 C[20] = fC[20];
2849 C[21] = fC[21] + dS*fC[24];
2850 C[22] = fC[22] + dS*fC[25];
2851 C[23] = fC[23] + dS*fC[26];
2852
2853 C[24] = fC[24];
2854 C[25] = fC[25];
2855 C[26] = fC[26];
2856 C[27] = fC[27];
2857 C[28] = fC[28] + dS*fC[31];
2858 C[29] = fC[29] + dS*fC[32];
2859 C[30] = fC[30] + dS*fC[33];
2860
2861 C[31] = fC[31];
2862 C[32] = fC[32];
2863 C[33] = fC[33];
2864 C[34] = fC[34];
2865 C[35] = fC[35];
2866}
2867
2869 const KFParticleBaseSIMD &daughter2, fvec Bz )
2870{
2871 //* Create gamma
2872
2873 const KFParticleBaseSIMD *daughters[2] = { &daughter1, &daughter2};
2874
2875 fvec v0[3];
2876
2877 if( !fIsLinearized ){
2878 fvec ds, ds1;
2879 fvec m[8];
2880 fvec mCd[36];
2881 daughter1.GetDStoParticle(daughter2, ds, ds1);
2882 daughter1.Transport( ds, m, mCd );
2883 fP[0] = m[0];
2884 fP[1] = m[1];
2885 fP[2] = m[2];
2886 daughter2.Transport( ds1, m, mCd );
2887 fP[0] = .5*( fP[0] + m[0] );
2888 fP[1] = .5*( fP[1] + m[1] );
2889 fP[2] = .5*( fP[2] + m[2] );
2890 } else {
2891 fP[0] = fVtxGuess[0];
2892 fP[1] = fVtxGuess[1];
2893 fP[2] = fVtxGuess[2];
2894 }
2895
2896 fvec daughterP[2][8], daughterC[2][36];
2897 fvec vtxMom[2][3];
2898
2899 int nIter = fIsLinearized ?1 :2;
2900
2901 for( int iter=0; iter<nIter; iter++){
2902
2903 v0[0] = fP[0];
2904 v0[1] = fP[1];
2905 v0[2] = fP[2];
2906
2908 fSFromDecay = 0;
2909 fP[0] = v0[0];
2910 fP[1] = v0[1];
2911 fP[2] = v0[2];
2912 fP[3] = 0;
2913 fP[4] = 0;
2914 fP[5] = 0;
2915 fP[6] = 0;
2916 fP[7] = 0;
2917
2918
2919 // fit daughters to the vertex guess
2920
2921 {
2922 for( int id=0; id<2; id++ ){
2923
2924 fvec *p = daughterP[id];
2925 fvec *mC = daughterC[id];
2926
2927 daughters[id]->GetMeasurement( v0, p, mC );
2928
2929 fvec mAi[6];
2930 InvertSym3(mC, mAi );
2931
2932 fvec mB[3][3];
2933
2934 mB[0][0] = mC[ 6]*mAi[0] + mC[ 7]*mAi[1] + mC[ 8]*mAi[3];
2935 mB[0][1] = mC[ 6]*mAi[1] + mC[ 7]*mAi[2] + mC[ 8]*mAi[4];
2936 mB[0][2] = mC[ 6]*mAi[3] + mC[ 7]*mAi[4] + mC[ 8]*mAi[5];
2937
2938 mB[1][0] = mC[10]*mAi[0] + mC[11]*mAi[1] + mC[12]*mAi[3];
2939 mB[1][1] = mC[10]*mAi[1] + mC[11]*mAi[2] + mC[12]*mAi[4];
2940 mB[1][2] = mC[10]*mAi[3] + mC[11]*mAi[4] + mC[12]*mAi[5];
2941
2942 mB[2][0] = mC[15]*mAi[0] + mC[16]*mAi[1] + mC[17]*mAi[3];
2943 mB[2][1] = mC[15]*mAi[1] + mC[16]*mAi[2] + mC[17]*mAi[4];
2944 mB[2][2] = mC[15]*mAi[3] + mC[16]*mAi[4] + mC[17]*mAi[5];
2945
2946 fvec z[3] = { v0[0]-p[0], v0[1]-p[1], v0[2]-p[2] };
2947
2948 vtxMom[id][0] = p[3] + mB[0][0]*z[0] + mB[0][1]*z[1] + mB[0][2]*z[2];
2949 vtxMom[id][1] = p[4] + mB[1][0]*z[0] + mB[1][1]*z[1] + mB[1][2]*z[2];
2950 vtxMom[id][2] = p[5] + mB[2][0]*z[0] + mB[2][1]*z[1] + mB[2][2]*z[2];
2951
2952 daughters[id]->Transport( daughters[id]->GetDStoPoint(v0), p, mC );
2953
2954 }
2955
2956 } // fit daughters to guess
2957
2958
2959 // fit new vertex
2960 {
2961
2962 fvec mpx0 = vtxMom[0][0]+vtxMom[1][0];
2963 fvec mpy0 = vtxMom[0][1]+vtxMom[1][1];
2964 fvec mpt0 = sqrt(mpx0*mpx0 + mpy0*mpy0);
2965 // fvec a0 = TMath::ATan2(mpy0,mpx0);
2966
2967 fvec ca0 = mpx0/mpt0;
2968 fvec sa0 = mpy0/mpt0;
2969 fvec r[3] = { v0[0], v0[1], v0[2] };
2970 fvec mC[3][3] = {{1000., 0 , 0 },
2971 {0, 1000., 0 },
2972 {0, 0, 1000. } };
2973 fvec chi2=0;
2974
2975 for( int id=0; id<2; id++ ){
2976 const fvec kCLight = 0.000299792458;
2977 fvec q = Bz*daughters[id]->GetQ()*kCLight;
2978 fvec px0 = vtxMom[id][0];
2979 fvec py0 = vtxMom[id][1];
2980 fvec pz0 = vtxMom[id][2];
2981 fvec pt0 = sqrt(px0*px0+py0*py0);
2982 fvec mG[3][6], mB[3], mH[3][3];
2983 // r = {vx,vy,vz};
2984 // m = {x,y,z,Px,Py,Pz};
2985 // V = daughter.C
2986 // G*m + B = H*r;
2987 // q*x + Py - q*vx - sin(a)*Pt = 0
2988 // q*y - Px - q*vy + cos(a)*Pt = 0
2989 // (Px*cos(a) + Py*sin(a) ) (vz -z) - Pz( cos(a)*(vx-x) + sin(a)*(vy-y)) = 0
2990
2991 mG[0][0] = q;
2992 mG[0][1] = 0;
2993 mG[0][2] = 0;
2994 mG[0][3] = -sa0*px0/pt0;
2995 mG[0][4] = 1 -sa0*py0/pt0;
2996 mG[0][5] = 0;
2997 mH[0][0] = q;
2998 mH[0][1] = 0;
2999 mH[0][2] = 0;
3000 mB[0] = py0 - sa0*pt0 - mG[0][3]*px0 - mG[0][4]*py0 ;
3001
3002 // q*y - Px - q*vy + cos(a)*Pt = 0
3003
3004 mG[1][0] = 0;
3005 mG[1][1] = q;
3006 mG[1][2] = 0;
3007 mG[1][3] = -1 + ca0*px0/pt0;
3008 mG[1][4] = + ca0*py0/pt0;
3009 mG[1][5] = 0;
3010 mH[1][0] = 0;
3011 mH[1][1] = q;
3012 mH[1][2] = 0;
3013 mB[1] = -px0 + ca0*pt0 - mG[1][3]*px0 - mG[1][4]*py0 ;
3014
3015 // (Px*cos(a) + Py*sin(a) ) (z -vz) - Pz( cos(a)*(x-vx) + sin(a)*(y-vy)) = 0
3016
3017 mG[2][0] = -pz0*ca0;
3018 mG[2][1] = -pz0*sa0;
3019 mG[2][2] = px0*ca0 + py0*sa0;
3020 mG[2][3] = 0;
3021 mG[2][4] = 0;
3022 mG[2][5] = 0;
3023
3024 mH[2][0] = mG[2][0];
3025 mH[2][1] = mG[2][1];
3026 mH[2][2] = mG[2][2];
3027
3028 mB[2] = 0;
3029
3030 // fit the vertex
3031
3032 // V = GVGt
3033
3034 fvec mGV[3][6];
3035 fvec mV[6];
3036 fvec m[3];
3037 for( int i=0; i<3; i++ ){
3038 m[i] = mB[i];
3039 for( int k=0; k<6; k++ ) m[i]+=mG[i][k]*daughterP[id][k];
3040 }
3041 for( int i=0; i<3; i++ ){
3042 for( int j=0; j<6; j++ ){
3043 mGV[i][j] = 0;
3044 for( int k=0; k<6; k++ ) mGV[i][j]+=mG[i][k]*daughterC[id][ IJ(k,j) ];
3045 }
3046 }
3047 for( int i=0, k=0; i<3; i++ ){
3048 for( int j=0; j<=i; j++,k++ ){
3049 mV[k] = 0;
3050 for( int l=0; l<6; l++ ) mV[k]+=mGV[i][l]*mG[j][l];
3051 }
3052 }
3053
3054
3055 //* CHt
3056
3057 fvec mCHt[3][3];
3058 fvec mHCHt[6];
3059 fvec mHr[3];
3060 for( int i=0; i<3; i++ ){
3061 mHr[i] = 0;
3062 for( int k=0; k<3; k++ ) mHr[i]+= mH[i][k]*r[k];
3063 }
3064
3065 for( int i=0; i<3; i++ ){
3066 for( int j=0; j<3; j++){
3067 mCHt[i][j] = 0;
3068 for( int k=0; k<3; k++ ) mCHt[i][j]+= mC[i][k]*mH[j][k];
3069 }
3070 }
3071
3072 for( int i=0, k=0; i<3; i++ ){
3073 for( int j=0; j<=i; j++, k++ ){
3074 mHCHt[k] = 0;
3075 for( int l=0; l<3; l++ ) mHCHt[k]+= mH[i][l]*mCHt[l][j];
3076 }
3077 }
3078
3079 fvec mS[6] = { mHCHt[0]+mV[0],
3080 mHCHt[1]+mV[1], mHCHt[2]+mV[2],
3081 mHCHt[3]+mV[3], mHCHt[4]+mV[4], mHCHt[5]+mV[5] };
3082
3083
3084 InvertSym3(mS,mS);
3085
3086 //* Residual (measured - estimated)
3087
3088 fvec zeta[3] = { m[0]-mHr[0], m[1]-mHr[1], m[2]-mHr[2] };
3089
3090 //* Kalman gain K = mCH'*S
3091
3092 fvec k[3][3];
3093
3094 for(Int_t i=0;i<3;++i){
3095 k[i][0] = mCHt[i][0]*mS[0] + mCHt[i][1]*mS[1] + mCHt[i][2]*mS[3];
3096 k[i][1] = mCHt[i][0]*mS[1] + mCHt[i][1]*mS[2] + mCHt[i][2]*mS[4];
3097 k[i][2] = mCHt[i][0]*mS[3] + mCHt[i][1]*mS[4] + mCHt[i][2]*mS[5];
3098 }
3099
3100 //* New estimation of the vertex position r += K*zeta
3101
3102 for(Int_t i=0;i<3;++i)
3103 r[i] = r[i] + k[i][0]*zeta[0] + k[i][1]*zeta[1] + k[i][2]*zeta[2];
3104
3105 //* New covariance matrix C -= K*(mCH')'
3106
3107 for(Int_t i=0;i<3;++i){
3108 for(Int_t j=0;j<=i;++j){
3109 mC[i][j] = mC[i][j] - (k[i][0]*mCHt[j][0] + k[i][1]*mCHt[j][1] + k[i][2]*mCHt[j][2]);
3110 mC[j][i] = mC[i][j];
3111 }
3112 }
3113
3114 //* Calculate Chi^2
3115
3116 chi2 += ( ( mS[0]*zeta[0] + mS[1]*zeta[1] + mS[3]*zeta[2] )*zeta[0]
3117 +(mS[1]*zeta[0] + mS[2]*zeta[1] + mS[4]*zeta[2] )*zeta[1]
3118 +(mS[3]*zeta[0] + mS[4]*zeta[1] + mS[5]*zeta[2] )*zeta[2] );
3119 }
3120
3121 // store vertex
3122
3123 fNDF = 2;
3124 fChi2 = chi2;
3125 for( int i=0; i<3; i++ ) fP[i] = r[i];
3126 for( int i=0,k=0; i<3; i++ ){
3127 for( int j=0; j<=i; j++,k++ ){
3128 fC[k] = mC[i][j];
3129 }
3130 }
3131 }
3132
3133 } // iterations
3134
3135 // now fit daughters to the vertex
3136
3137 fQ = 0;
3138 fSFromDecay = 0;
3139
3140 for(Int_t i=3;i<8;++i) fP[i]=0.;
3141 for(Int_t i=6;i<35;++i) fC[i]=0.;
3142 fC[35] = 100.;
3143
3144 for( int id=0; id<2; id++ ){
3145
3146 fvec *p = daughterP[id];
3147 fvec *mC = daughterC[id];
3148 daughters[id]->GetMeasurement( v0, p, mC );
3149
3150 const fvec *m = fP, *mV = fC;
3151
3152 fvec mAi[6];
3153 InvertSym3(mC, mAi );
3154
3155 fvec mB[4][3];
3156
3157 mB[0][0] = mC[ 6]*mAi[0] + mC[ 7]*mAi[1] + mC[ 8]*mAi[3];
3158 mB[0][1] = mC[ 6]*mAi[1] + mC[ 7]*mAi[2] + mC[ 8]*mAi[4];
3159 mB[0][2] = mC[ 6]*mAi[3] + mC[ 7]*mAi[4] + mC[ 8]*mAi[5];
3160
3161 mB[1][0] = mC[10]*mAi[0] + mC[11]*mAi[1] + mC[12]*mAi[3];
3162 mB[1][1] = mC[10]*mAi[1] + mC[11]*mAi[2] + mC[12]*mAi[4];
3163 mB[1][2] = mC[10]*mAi[3] + mC[11]*mAi[4] + mC[12]*mAi[5];
3164
3165 mB[2][0] = mC[15]*mAi[0] + mC[16]*mAi[1] + mC[17]*mAi[3];
3166 mB[2][1] = mC[15]*mAi[1] + mC[16]*mAi[2] + mC[17]*mAi[4];
3167 mB[2][2] = mC[15]*mAi[3] + mC[16]*mAi[4] + mC[17]*mAi[5];
3168
3169 mB[3][0] = mC[21]*mAi[0] + mC[22]*mAi[1] + mC[23]*mAi[3];
3170 mB[3][1] = mC[21]*mAi[1] + mC[22]*mAi[2] + mC[23]*mAi[4];
3171 mB[3][2] = mC[21]*mAi[3] + mC[22]*mAi[4] + mC[23]*mAi[5];
3172
3173
3174 fvec z[3] = { m[0]-p[0], m[1]-p[1], m[2]-p[2] };
3175
3176// {
3177// fvec mAV[6] = { mC[0]-mV[0], mC[1]-mV[1], mC[2]-mV[2],
3178// mC[3]-mV[3], mC[4]-mV[4], mC[5]-mV[5] };
3179//
3180// fvec mAVi[6];
3181// if( !InvertSym3(mAV, mAVi) ){
3182// fvec dChi2 = ( +(mAVi[0]*z[0] + mAVi[1]*z[1] + mAVi[3]*z[2])*z[0]
3183// +(mAVi[1]*z[0] + mAVi[2]*z[1] + mAVi[4]*z[2])*z[1]
3184// +(mAVi[3]*z[0] + mAVi[4]*z[1] + mAVi[5]*z[2])*z[2] );
3185// fChi2+= TMath::Abs( dChi2 );
3186// }
3187// fNDF += 2;
3188// }
3189
3190 //* Add the daughter momentum to the particle momentum
3191
3192 fP[3]+= p[3] + mB[0][0]*z[0] + mB[0][1]*z[1] + mB[0][2]*z[2];
3193 fP[4]+= p[4] + mB[1][0]*z[0] + mB[1][1]*z[1] + mB[1][2]*z[2];
3194 fP[5]+= p[5] + mB[2][0]*z[0] + mB[2][1]*z[1] + mB[2][2]*z[2];
3195 fP[6]+= p[6] + mB[3][0]*z[0] + mB[3][1]*z[1] + mB[3][2]*z[2];
3196
3197 fvec d0, d1, d2;
3198
3199 d0= mB[0][0]*mV[0] + mB[0][1]*mV[1] + mB[0][2]*mV[3] - mC[ 6];
3200 d1= mB[0][0]*mV[1] + mB[0][1]*mV[2] + mB[0][2]*mV[4] - mC[ 7];
3201 d2= mB[0][0]*mV[3] + mB[0][1]*mV[4] + mB[0][2]*mV[5] - mC[ 8];
3202
3203 //fC[6]+= mC[ 6] + d0;
3204 //fC[7]+= mC[ 7] + d1;
3205 //fC[8]+= mC[ 8] + d2;
3206 fC[9]+= mC[ 9] + d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
3207
3208 d0= mB[1][0]*mV[0] + mB[1][1]*mV[1] + mB[1][2]*mV[3] - mC[10];
3209 d1= mB[1][0]*mV[1] + mB[1][1]*mV[2] + mB[1][2]*mV[4] - mC[11];
3210 d2= mB[1][0]*mV[3] + mB[1][1]*mV[4] + mB[1][2]*mV[5] - mC[12];
3211
3212 //fC[10]+= mC[10]+ d0;
3213 //fC[11]+= mC[11]+ d1;
3214 //fC[12]+= mC[12]+ d2;
3215 fC[13]+= mC[13]+ d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
3216 fC[14]+= mC[14]+ d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
3217
3218 d0= mB[2][0]*mV[0] + mB[2][1]*mV[1] + mB[2][2]*mV[3] - mC[15];
3219 d1= mB[2][0]*mV[1] + mB[2][1]*mV[2] + mB[2][2]*mV[4] - mC[16];
3220 d2= mB[2][0]*mV[3] + mB[2][1]*mV[4] + mB[2][2]*mV[5] - mC[17];
3221
3222 //fC[15]+= mC[15]+ d0;
3223 //fC[16]+= mC[16]+ d1;
3224 //fC[17]+= mC[17]+ d2;
3225 fC[18]+= mC[18]+ d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
3226 fC[19]+= mC[19]+ d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
3227 fC[20]+= mC[20]+ d0*mB[2][0] + d1*mB[2][1] + d2*mB[2][2];
3228
3229 d0= mB[3][0]*mV[0] + mB[3][1]*mV[1] + mB[3][2]*mV[3] - mC[21];
3230 d1= mB[3][0]*mV[1] + mB[3][1]*mV[2] + mB[3][2]*mV[4] - mC[22];
3231 d2= mB[3][0]*mV[3] + mB[3][1]*mV[4] + mB[3][2]*mV[5] - mC[23];
3232
3233 //fC[21]+= mC[21] + d0;
3234 //fC[22]+= mC[22] + d1;
3235 //fC[23]+= mC[23] + d2;
3236 fC[24]+= mC[24] + d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
3237 fC[25]+= mC[25] + d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
3238 fC[26]+= mC[26] + d0*mB[2][0] + d1*mB[2][1] + d2*mB[2][2];
3239 fC[27]+= mC[27] + d0*mB[3][0] + d1*mB[3][1] + d2*mB[3][2];
3240 }
3241
3242// SetMassConstraint(0,0);
3244}
3245
3247{
3248// example:
3249// KFParticle PosParticle(...)
3250// KFParticle NegParticle(...)
3251// Gamma.ConstructGamma(PosParticle, NegParticle);
3252// fvec VertexGamma[3] = {Gamma.GetX(), Gamma.GetY(), Gamma.GetZ()};
3253// PosParticle.TransportToPoint(VertexGamma);
3254// NegParticle.TransportToPoint(VertexGamma);
3255// fvec armenterosQtAlfa[2] = {0.};
3256// KFParticle::GetArmenterosPodolanski(PosParticle, NegParticle, armenterosQtAlfa );
3257
3258 fvec alpha = 0., qt = 0.;
3259 fvec spx = positive.GetPx() + negative.GetPx();
3260 fvec spy = positive.GetPy() + negative.GetPy();
3261 fvec spz = positive.GetPz() + negative.GetPz();
3262 fvec sp = sqrt(spx*spx + spy*spy + spz*spz);
3263 fvec mask = fvec( fabs(sp) < fvec(1.E-10));
3264 fvec pn, pp, pln, plp;
3265
3266 pn = sqrt(negative.GetPx()*negative.GetPx() + negative.GetPy()*negative.GetPy() + negative.GetPz()*negative.GetPz());
3267 pp = sqrt(positive.GetPx()*positive.GetPx() + positive.GetPy()*positive.GetPy() + positive.GetPz()*positive.GetPz());
3268 pln = (negative.GetPx()*spx+negative.GetPy()*spy+negative.GetPz()*spz)/sp;
3269 plp = (positive.GetPx()*spx+positive.GetPy()*spy+positive.GetPz()*spz)/sp;
3270
3271 mask = fvec(mask & fvec( fabs(pn) < fvec(1.E-10)));
3272 fvec ptm = (1.-((pln/pn)*(pln/pn)));
3273 qt= if3( fvec(ptm >= Zero) , pn*sqrt(ptm) , Zero );
3274 alpha = (plp-pln)/(plp+pln);
3275
3276 QtAlfa[0] = fvec(qt & mask);
3277 QtAlfa[1] = fvec(alpha & mask);
3278}
3279
3281{
3282 // Rotates the KFParticle object around OZ axis, OZ axis is set by the vertex position
3283 // fvec angle - angle of rotation in XY plane in [rad]
3284 // fvec Vtx[3] - position of the vertex in [cm]
3285
3286 // Before rotation the center of the coordinat system should be moved to the vertex position; move back after rotation
3287 X() = X() - Vtx[0];
3288 Y() = Y() - Vtx[1];
3289 Z() = Z() - Vtx[2];
3290
3291 // Rotate the kf particle
3292 fvec c = cos(angle);
3293 fvec s = sin(angle);
3294
3295 fvec mA[8][ 8];
3296 for( Int_t i=0; i<8; i++ ){
3297 for( Int_t j=0; j<8; j++){
3298 mA[i][j] = 0;
3299 }
3300 }
3301 for( int i=0; i<8; i++ ){
3302 mA[i][i] = 1;
3303 }
3304 mA[0][0] = c; mA[0][1] = s;
3305 mA[1][0] = -s; mA[1][1] = c;
3306 mA[3][3] = c; mA[3][4] = s;
3307 mA[4][3] = -s; mA[4][4] = c;
3308
3309 fvec mAC[8][8];
3310 fvec mAp[8];
3311
3312 for( Int_t i=0; i<8; i++ ){
3313 mAp[i] = 0;
3314 for( Int_t k=0; k<8; k++){
3315 mAp[i]+=mA[i][k] * fP[k];
3316 }
3317 }
3318
3319 for( Int_t i=0; i<8; i++){
3320 fP[i] = mAp[i];
3321 }
3322
3323 for( Int_t i=0; i<8; i++ ){
3324 for( Int_t j=0; j<8; j++ ){
3325 mAC[i][j] = 0;
3326 for( Int_t k=0; k<8; k++ ){
3327 mAC[i][j]+= mA[i][k] * GetCovariance(k,j);
3328 }
3329 }
3330 }
3331
3332 for( Int_t i=0; i<8; i++ ){
3333 for( Int_t j=0; j<=i; j++ ){
3334 fvec xx = 0;
3335 for( Int_t k=0; k<8; k++){
3336 xx+= mAC[i][k]*mA[j][k];
3337 }
3338 Covariance(i,j) = xx;
3339 }
3340 }
3341
3342 X() = GetX() + Vtx[0];
3343 Y() = GetY() + Vtx[1];
3344 Z() = GetZ() + Vtx[2];
3345}
3346
3348{
3349 //* Invert symmetric matric stored in low-triagonal form
3350
3351 fvec ret = 0;
3352 fvec a0 = A[0], a1 = A[1], a2 = A[2], a3 = A[3];
3353
3354 Ai[0] = a2*A[5] - A[4]*A[4];
3355 Ai[1] = a3*A[4] - a1*A[5];
3356 Ai[3] = a1*A[4] - a2*a3;
3357 fvec det = (a0*Ai[0] + a1*Ai[1] + a3*Ai[3]);
3358
3359 fvec idet = 1. / det;
3360 fvec init = fvec(small < fabs(det));
3361 det = init & idet;
3362 ret += (!init);
3363
3364 Ai[0] *= det;
3365 Ai[1] *= det;
3366 Ai[3] *= det;
3367 Ai[2] = ( a0*A[5] - a3*a3 )*det;
3368 Ai[4] = ( a1*a3 - a0*A[4] )*det;
3369 Ai[5] = ( a0*a2 - a1*a1 )*det;
3370 return ret;
3371}
3372
3374{
3375 fvec d[3], uud, u[3][3];
3376 for(int i=0; i<3; i++)
3377 {
3378 d[i]=Zero;
3379 for(int j=0; j<3; j++)
3380 u[i][j]=Zero;
3381 }
3382
3383 for(int i=0; i<3 ; i++)
3384 {
3385 uud=Zero;
3386 for(int j=0; j<i; j++)
3387 uud += u[j][i]*u[j][i]*d[j];
3388 uud = a[i*(i+3)/2] - uud;
3389
3390 fvec smallval = 1.e-12f;
3391 fvec initialised = fvec( fabs(uud)<smallval );
3392 uud = ((!initialised) & uud) + (smallval & initialised);
3393
3394 d[i] = uud/fabs(uud);
3395 u[i][i] = sqrt(fabs(uud));
3396
3397 for(int j=i+1; j<3; j++)
3398 {
3399 uud = Zero;
3400 for(int k=0; k<i; k++)
3401 uud += u[k][i]*u[k][j]*d[k];
3402 uud = a[j*(j+1)/2+i] - uud;
3403 u[i][j] = d[i]/u[i][i]*uud;
3404 }
3405 }
3406
3407 fvec u1[3];
3408
3409 for(int i=0; i<3; i++)
3410 {
3411 u1[i] = u[i][i];
3412 u[i][i] = One/u[i][i];
3413 }
3414 for(int i=0; i<2; i++)
3415 {
3416 u[i][i+1] = - u[i][i+1]*u[i][i]*u[i+1][i+1];
3417 }
3418 for(int i=0; i<1; i++)
3419 {
3420 u[i][i+2] = u[i][i+1]*u1[i+1]*u[i+1][i+2]-u[i][i+2]*u[i][i]*u[i+2][i+2];
3421 }
3422
3423 for(int i=0; i<3; i++)
3424 a[i+3] = u[i][2]*u[2][2]*d[2];
3425 for(int i=0; i<2; i++)
3426 a[i+1] = u[i][1]*u[1][1]*d[1] + u[i][2]*u[1][2]*d[2];
3427 a[0] = u[0][0]*u[0][0]*d[0] + u[0][1]*u[0][1]*d[1] + u[0][2]*u[0][2]*d[2];
3428}
3429
3430void KFParticleBaseSIMD::MultQSQt( const fvec Q[], const fvec S[], fvec SOut[] )
3431{
3432 //* Matrix multiplication Q*S*Q^T, Q - square matrix, S - symmetric
3433
3434 const Int_t kN= 8;
3435 fvec mA[kN*kN];
3436
3437 for( Int_t i=0, ij=0; i<kN; i++ ){
3438 for( Int_t j=0; j<kN; j++, ++ij ){
3439 mA[ij] = 0 ;
3440 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];
3441 }
3442 }
3443
3444 for( Int_t i=0; i<kN; i++ ){
3445 for( Int_t j=0; j<=i; j++ ){
3446 Int_t ij = ( j<=i ) ? i*(i+1)/2+j :j*(j+1)/2+i;
3447 SOut[ij] = 0 ;
3448 for( Int_t k=0; k<kN; k++ ) SOut[ij] += Q[ i*kN+k ] * mA[ k*kN+j ];
3449 }
3450 }
3451}
3452
3453
3454// 72-charachters line to define the printer border
3455//3456789012345678901234567890123456789012345678901234567890123456789012
3456
const Float_t d
Z-ccordinate of the first GEM-station.
Definition BmnMwpcHit.cxx:7
friend F32vec4 sin(const F32vec4 &a)
Definition P4_F32vec4.h:124
friend F32vec4 sqrt(const F32vec4 &a)
Definition P4_F32vec4.h:34
friend F32vec4 fabs(const F32vec4 &a)
Definition P4_F32vec4.h:52
__m128 v
Definition P4_F32vec4.h:1
friend F32vec4 cos(const F32vec4 &a)
Definition P4_F32vec4.h:125
int i
Definition P4_F32vec4.h:22
__m128 m
Definition P4_F32vec4.h:27
friend F32vec4 atan2(const F32vec4 &y, const F32vec4 &x)
Definition P4_F32vec4.h:130
#define if3(a, b, c)
Definition P4_F32vec4.h:148
friend F32vec4 log(const F32vec4 &a)
Definition P4_F32vec4.h:123
float f
Definition P4_F32vec4.h:21
F32vec4 fvec
Definition P4_F32vec4.h:231
void GetDStoParticleBz(fvec Bz, const KFParticleBaseSIMD &p, fvec &dS, fvec &dS1) const
static void InvertCholetsky3(fvec a[6])
static fvec InvertSym3(const fvec A[], fvec Ainv[])
void AddDaughterWithEnergyFitMC(const KFParticleBaseSIMD &Daughter, Bool_t isAtVtxGuess)
virtual void GetFieldValue(const fvec xyz[], fvec B[]) const =0
void SubtractFromParticle(KFParticleBaseSIMD &Vtx) const
const fvec & Q() const
fvec & Covariance(Int_t i)
void ConstructGammaBz(const KFParticleBaseSIMD &daughter1, const KFParticleBaseSIMD &daughter2, fvec Bz)
void TransportLine(fvec S, fvec P[], fvec C[]) const
fvec GetDeviationFromParticle(const KFParticleBaseSIMD &p) const
void GetDistanceToVertexLine(const KFParticleBaseSIMD &Vertex, fvec &l, fvec &dl, fvec *isParticleFromVertex=0) const
static Int_t IJ(Int_t i, Int_t j)
void GetMeasurement(const fvec XYZ[], fvec m[], fvec V[], Bool_t isAtVtxGuess=0) const
fvec GetDistanceFromParticle(const KFParticleBaseSIMD &p) const
void operator+=(const KFParticleBaseSIMD &Daughter)
fvec GetDStoPointBy(fvec By, const fvec xyz[]) const
fvec GetMomentum(fvec &P, fvec &SigmaP) const
const fvec & S() const
void Convert(bool ToProduction)
void SetMassConstraint(fvec Mass, fvec SigmaMass=0)
fvec GetR(fvec &R, fvec &SigmaR) const
void AddDaughter(const KFParticleBaseSIMD &Daughter, Bool_t isAtVtxGuess=0)
void AddDaughterWithEnergyFit(const KFParticleBaseSIMD &Daughter, Bool_t isAtVtxGuess)
virtual fvec GetDStoPoint(const fvec xyz[]) const =0
fvec GetSCorrection(const fvec Part[], const fvec XYZ[]) const
fvec GetDistanceFromVertex(const fvec vtx[]) const
const fvec & X() const
void TransportCBM(fvec dS, fvec P[], fvec C[]) const
fvec GetDecayLength(fvec &L, fvec &SigmaL) const
void GetDStoParticleLine(const KFParticleBaseSIMD &p, fvec &dS, fvec &dS1) const
fvec GetDecayLengthXY(fvec &L, fvec &SigmaL) const
void RotateXY(fvec angle, fvec Vtx[3])
void SetVtxGuess(fvec x, fvec y, fvec z)
static void multQSQt1(const fvec J[11], fvec S[])
void SetNonlinearMassConstraint(fvec Mass)
virtual void GetDStoParticle(const KFParticleBaseSIMD &p, fvec &DS, fvec &DSp) const =0
fvec GetCovariance(Int_t i) const
fvec GetDeviationFromVertex(const fvec v[], const fvec Cv[]=0) const
virtual void Transport(fvec dS, fvec P[], fvec C[]) const =0
fvec GetLifeTime(fvec &T, fvec &SigmaT) const
void GetDStoParticleCBM(const KFParticleBaseSIMD &p, fvec &dS, fvec &dS1) const
fvec GetDStoPointCBM(const fvec xyz[]) const
fvec GetDStoPointBz(fvec Bz, const fvec xyz[]) const
const fvec & Y() const
void SetVtxErrGuess(fvec &x, fvec &y, fvec &z)
void GetDStoParticleBy(fvec B, const KFParticleBaseSIMD &p, fvec &dS, fvec &dS1) const
void SetProductionVertex(const KFParticleBaseSIMD &Vtx)
void Construct(const KFParticleBaseSIMD *vDaughters[], Int_t nDaughters, const KFParticleBaseSIMD *ProdVtx=0, Float_t Mass=-1, Bool_t IsConstrained=0, Bool_t isAtVtxGuess=0)
void AddDaughterWithEnergyCalc(const KFParticleBaseSIMD &Daughter, Bool_t isAtVtxGuess)
fvec GetPt(fvec &Pt, fvec &SigmaPt) const
void SubtractFromVertex(KFParticleBaseSIMD &Vtx) const
void TransportBz(fvec Bz, fvec dS, fvec P[], fvec C[]) const
fvec GetPhi(fvec &Phi, fvec &SigmaPhi) const
fvec & Cij(Int_t i, Int_t j)
fvec GetEta(fvec &Eta, fvec &SigmaEta) const
fvec GetMass(fvec &M, fvec &SigmaM) const
const fvec & Z() const
static void MultQSQt(const fvec Q[], const fvec S[], fvec SOut[])
static void GetArmenterosPodolanski(KFParticleBaseSIMD &positive, KFParticleBaseSIMD &negative, fvec QtAlfa[2])
void AddDaughterId(fvec id)