BmnRoot
Loading...
Searching...
No Matches
CbmKFSecondaryVertexFinder.cxx
Go to the documentation of this file.
1
10
11#include "CbmKF.h"
12#include "CbmKFTrack.h"
13
14#include <cmath>
15
16using std::vector;
17using std::fabs;
18
20{
21 vTracks.clear();
22 VGuess =0;
23 VParent = 0;
24 MassConstraint = -1;
25}
26
28{
29 vTracks.clear();
30}
31
33{
34 vTracks.push_back(Track);
35}
36
37void CbmKFSecondaryVertexFinder::SetTracks( vector<CbmKFTrackInterface*> &vTr )
38{
39 vTracks = vTr;
40}
41
46
48{
49 MassConstraint = MotherMass;
50}
51
56
58{
59 const Int_t MaxIter=3;
60
61 if( VGuess ){
62 r[0] = VGuess->GetRefX();
63 r[1] = VGuess->GetRefY();
64 r[2] = VGuess->GetRefZ();
65 }else {
66 if( CbmKF::Instance()->vTargets.empty() ){
67 r[0] = r[1] = r[2] = 0.;
68 }else{
70 r[0] = r[1] = 0.;
71 r[2] = t.z;
72 }
73 }
74
75 r[3] = 0;
76 r[4] = 0;
77 r[5] = 0;
78 r[6] = 0;
79 r[7] = 0;
80
81 for ( Int_t iteration =0; iteration<MaxIter;++iteration ){
82
83 for( int i=0; i<8; i++) r0[i] = r[i];
84
85 r[3] = 0;
86 r[4] = 0;
87 r[5] = 0;
88 r[6] = 0;
89
90 for(Int_t i=0;i<36;++i) C[i]=0.0;
91
92 C[0] = C[2] = C[5] = 100.0;
93 C[35] = 1.E4;
94 NDF = -3;
95 Chi2 = 0.;
96
97 for( vector<CbmKFTrackInterface*>::iterator tr=vTracks.begin(); tr!=vTracks.end(); ++tr )
98 {
99 CbmKFTrack T( **tr );
100
101 T.Extrapolate( r0[2] );
102
103 Double_t *m = T.GetTrack();
104 Double_t *V = T.GetCovMatrix();
105
106 //* Fit of vertex track slopes and momentum (a,b,qp) to r0 vertex
107
108 Double_t a = m[2], b = m[3], qp = m[4];
109 {
110 Double_t s = V[0]*V[2] - V[1]*V[1];
111 s = ( s > 1.E-20 ) ?1./s :0;
112 Double_t zetas[2] = { (r0[0]-m[0])*s, (r0[1]-m[1])*s };
113 a += ( V[ 3]*V[2] -V[ 4]*V[1])*zetas[0]
114 +(-V[ 3]*V[1] +V[ 4]*V[0])*zetas[1];
115 b += ( V[ 6]*V[2] -V[ 7]*V[1])*zetas[0]
116 +(-V[ 6]*V[1] +V[ 7]*V[0])*zetas[1];
117 qp+= (V[10]*V[2]-V[11]*V[1])*zetas[0]
118 +(-V[10]*V[1]+V[11]*V[0])*zetas[1];
119 }
120
121 //* convert the track to (x,y,px,py,pz,e) parameterization
122 double mm[6], VV[21];
123 {
124 double c2 = 1./(1. + a*a + b*b);
125 double p = 1./qp;
126 double p2 = p*p;
127 double pz = sqrt(p2*c2);
128 double px = a*pz;
129 double py = b*pz;
130 double E = sqrt( T.GetMass()*T.GetMass() + p2 );
131
132 double da = (m[2]-a);
133 double db = (m[3]-b);
134 double dq = (m[4]-qp);
135
136 double H[3] = { -px*c2, -py*c2, -pz*p };
137 double HE = -p*p2/E;
138
139 double dpz = H[0]*da + H[1]*db + H[2]*dq;
140
141 mm[0] = m[0];
142 mm[1] = m[1];
143 mm[2] = px + da + a*dpz;
144 mm[3] = py + db + b*dpz;
145 mm[4] = pz + dpz;
146 mm[5] = E + HE*dq;
147
148 double cxpz = H[0]*V[ 3] + H[1]*V[ 6] + H[2]*V[10];
149 double cypz = H[0]*V[ 4] + H[1]*V[ 7] + H[2]*V[11];
150 double capz = H[0]*V[ 5] + H[1]*V[ 8] + H[2]*V[12];
151 double cbpz = H[0]*V[ 8] + H[1]*V[ 9] + H[2]*V[13];
152 double cqpz = H[0]*V[12] + H[1]*V[13] + H[2]*V[14];
153 double cpzpz = H[0]*H[0]*V[5] +H[1]*H[1]*V[9] + H[2]*H[2]*V[14]
154 + 2*( H[0]*H[1]*V[8] +H[0]*H[2]*V[12] +H[1]*H[2]*V[13]);
155
156 VV[ 0] = V[0];
157 VV[ 1] = V[1];
158 VV[ 2] = V[2];
159 VV[ 3] = V[3]*pz + a*cxpz;
160 VV[ 4] = V[4]*pz + a*cypz;
161 VV[ 5] = V[5]*pz*pz + 2*a*pz*capz + a*a*cpzpz;
162 VV[ 6] = V[6]*pz+b*cxpz;
163 VV[ 7] = V[7]*pz+b*cypz;
164 VV[ 8] = V[8]*pz*pz + a*pz*cbpz + b*pz*capz + a*b*cpzpz;
165 VV[ 9] = V[9]*pz*pz + 2*b*pz*cbpz + b*b*cpzpz;
166 VV[10] = cxpz;
167 VV[11] = cypz;
168 VV[12] = capz*pz + a*cpzpz;
169 VV[13] = cbpz*pz + b*cpzpz;
170 VV[14] = cpzpz;
171 VV[15] = HE*V[10];
172 VV[16] = HE*V[11];
173 VV[17] = HE*( V[12]*pz + a*cqpz );
174 VV[18] = HE*( V[13]*pz + b*cqpz );
175 VV[19] = HE*cqpz;
176 VV[20] = HE*HE*V[14];
177 }
178
179
180 //* Measure the state vector with the track estimate
181
182 Chi2 += T.GetRefChi2();
183 NDF += T.GetRefNDF() - 3;
184
185 //* Update the state vector with the momentum part of the track estimate
186
187 r[ 3] += mm[ 2];
188 r[ 4] += mm[ 3];
189 r[ 5] += mm[ 4];
190 r[ 6] += mm[ 5];
191
192 C[ 9] += VV[ 5];
193 C[13] += VV[ 8];
194 C[14] += VV[ 9];
195 C[18] += VV[12];
196 C[19] += VV[13];
197 C[20] += VV[14];
198 C[24] += VV[17];
199 C[25] += VV[18];
200 C[26] += VV[19];
201 C[27] += VV[20];
202
203 NDF += 3;
204 Chi2 += 0.;
205
206 //* Update the state vector with the coordinate part of the track estimate
207
208 //* Residual (measured - estimated)
209
210 Double_t zeta[2] = { mm[0]-(r[0]-a*(r[2]-r0[2])),
211 mm[1]-(r[1]-b*(r[2]-r0[2])) };
212
213 //* Measurement matrix H= { { 1, 0, -a, 0..0}, { 0, 1, -b,0..0} };
214
215 //* S = (H*C*H' + V )^{-1}
216
217 Double_t S[3] = { VV[0] + C[0] - 2*a*C[3] + a*a*C[5],
218 VV[1] + C[1] - b*C[3] - a*C[4] + a*b*C[5],
219 VV[2] + C[2] - 2*b*C[4] + b*b*C[5] };
220
221 //* Invert S
222 {
223 Double_t s = S[0]*S[2] - S[1]*S[1];
224
225 if ( s < 1.E-20 ) continue;
226 s = 1./s;
227 Double_t S0 = S[0];
228 S[0] = s*S[2];
229 S[1] = -s*S[1];
230 S[2] = s*S0;
231 }
232
233 //* CHt = CH' - D'
234
235 Double_t CHt0[7], CHt1[7];
236
237 CHt0[0]= C[ 0] -a*C[ 3]; CHt1[0]= C[ 1] -b*C[ 3];
238 CHt0[1]= C[ 1] -a*C[ 4]; CHt1[1]= C[ 2] -b*C[ 4];
239 CHt0[2]= C[ 3] -a*C[ 5]; CHt1[2]= C[ 4] -b*C[ 5];
240 CHt0[3]= C[ 6] -a*C[ 8] -VV[ 3]; CHt1[3]= C[ 7] -b*C[ 8] -VV[ 4];
241 CHt0[4]= C[10] -a*C[12] -VV[ 6]; CHt1[4]= C[11] -b*C[12] -VV[ 7];
242 CHt0[5]= C[15] -a*C[17] -VV[10]; CHt1[5]= C[16] -b*C[17] -VV[11];
243 CHt0[6]= C[21] -a*C[23] -VV[15]; CHt1[6]= C[22] -b*C[23] -VV[16];
244
245 //* Kalman gain K = CH'*S
246
247 Double_t K0[7], K1[7];
248
249 for(Int_t i=0;i<7;++i){
250 K0[i] = CHt0[i]*S[0] + CHt1[i]*S[1];
251 K1[i] = CHt0[i]*S[1] + CHt1[i]*S[2];
252 }
253
254 //* New estimation of the vertex position r += K*zeta
255
256 for(Int_t i=0;i<7;++i) r[i] += K0[i]*zeta[0] + K1[i]*zeta[1];
257
258 //* New covariance matrix C -= K*(CH')'
259
260 for(Int_t i=0, k=0;i<7;++i){
261 for(Int_t j=0;j<=i;++j,++k)
262 C[k] -= K0[i]*CHt0[j] + K1[i]*CHt1[j];
263 }
264
265 //* Calculate Chi^2 & NDF
266
267 NDF += 2;
268 Chi2 += zeta[0]*zeta[0]*S[0] + 2*zeta[0]*zeta[1]*S[1]
269 + zeta[1]*zeta[1]*S[2] ;
270
271 } // add tracks
272
273 // Put constraints if they exist
274
275 AddTopoConstraint();
276 AddMassConstraint();
277
278 }// iterations
279}
280
281
283{
284 vtx.GetRefX() = r[0];
285 vtx.GetRefY() = r[1];
286 vtx.GetRefZ() = r[2];
287 for( int i=0; i<6; i++) vtx.GetCovMatrix()[i] = C[i];
288 vtx.GetRefChi2() = Chi2;
289 vtx.GetRefNDF() = NDF;
290}
291
293{
295 GetVertex( vi );
296 vi.GetVertex(vtx);
297}
298
299
301 GetMotherTrack( Double_t Par[], Double_t Cov[] )
302{
303
304 if( !Par ) return;
305
306 Double_t px=r[3], py=r[4], pz=r[5];
307
308 Double_t p = sqrt( px*px + py*py + pz*pz );
309 double pzi = 1/pz;
310 double qp=1/p;
311 Par[0] = r[0];
312 Par[1] = r[1];
313 Par[2] = px*pzi;
314 Par[3] = py*pzi;
315 Par[4] = qp;
316 Par[5] = r[2];
317
318 if( !Cov ) return;
319
320 double qp3=qp*qp*qp;
321 Double_t J[5][6];
322
323 for ( Int_t i=0; i<5; i++)
324 for (Int_t j=0;j<6;++j) J[i][j]=0.0;
325
326
327 J[0][0] = 1; J[0][2] = -Par[2];
328 J[1][1] = 1; J[1][2] = -Par[3];
329 J[2][3] = pzi; J[2][5] = -Par[2]*pzi;
330 J[3][4] = pzi; J[3][5] = -Par[3]*pzi;
331 J[4][3] = -qp3*px; J[4][4] = -qp3*py; J[4][4] = -qp3*pz;
332
333 Double_t JC[5][6];
334
335 for( Int_t i=0; i<5; i++)
336 for(Int_t j=0;j<6;j++){
337 JC[i][j]=0;
338 for( Int_t k=0; k<6; k++) JC[i][j]+=J[i][k]*Cij(k,j);
339 }
340
341 Int_t ii =0 ;
342 for( Int_t i=0; i<5; i++)
343 for(Int_t j=i;j<5;j++, ii++){
344 Cov[ii]=0;
345 for( Int_t k=0; k<6; k++) Cov[ii]+=JC[i][k]*J[j][k];
346 }
347}
348
349
350void CbmKFSecondaryVertexFinder::GetMass( Double_t *M, Double_t *Error_ )
351{
352 // S = sigma^2 of m2/2
353
354 Double_t S = ( r[3]*r[3]*C[9] + r[4]*r[4]*C[14] + r[5]*r[5]*C[20]
355 + r[6]*r[6]*C[27]
356 + 2*( + r[3]*r[4]*C[13] + r[5]*(r[3]*C[18] + r[4]*C[19])
357 - r[6]*( r[3]*C[24] + r[4]*C[25] + r[5]*C[26] ) )
358 );
359 Double_t m2 = r[6]*r[6] - r[3]*r[3] - r[4]*r[4] - r[5]*r[5];
360
361 if( M ) *M = ( m2>1.e-20 ) ? sqrt(m2) :0. ;
362 if( Error_ ) *Error_ = ( S>=0 && m2>1.e-20 ) ? sqrt(S/m2) :1.e4;
363}
364
365
366void CbmKFSecondaryVertexFinder::AddMassConstraint()
367{
368 if( MassConstraint<0 ) return;
369 double H[8];
370 H[0] = H[1] = H[2] = 0.;
371 H[3] = -2*r0[3];
372 H[4] = -2*r0[4];
373 H[5] = -2*r0[5];
374 H[6] = 2*r0[6];
375 H[7] = 0;
376 double m2 = r0[6]*r0[6] - r0[3]*r0[3] - r0[4]*r0[4] - r0[5]*r0[5];
377
378 double zeta = MassConstraint*MassConstraint - m2;
379 for(Int_t i=0;i<8;++i) zeta -= H[i]*(r[i]-r0[i]);
380
381 Double_t S = 0.;
382 Double_t CHt[8];
383 for (Int_t i=0;i<8;++i ){
384 CHt[i] = 0.0;
385 for (Int_t j=0;j<8;++j) CHt[i] += Cij(i,j)*H[j];
386 S += H[i]*CHt[i];
387 }
388
389 if( S<1.e-20 ) return;
390 S = 1./S;
391 Chi2 += zeta*zeta*S;
392 NDF += 1;
393 for( Int_t i=0, ii=0; i<8; ++i ){
394 Double_t Ki = CHt[i]*S;
395 r[i]+= Ki*zeta;
396 for(Int_t j=0;j<=i;++j) C[ii++] -= Ki*CHt[j];
397 }
398}
399
400
401
402
404{
405
406 r0[0]+= T*r0[3];
407 r0[1]+= T*r0[4];
408 r0[2]+= T*r0[5];
409
410 r[0]+= T*r[3];
411 r[1]+= T*r[4];
412 r[2]+= T*r[5];
413
414 double C30 = C[ 6] + T*C[ 9];
415 double C41 = C[11] + T*C[14];
416 double C52 = C[17] + T*C[20];
417 double T54 = T*C[19];
418
419 C[ 0] += T*( C[ 6] + C30 );
420 C[ 1] += T*( C[ 7] + C[10] );
421 C[ 2] += T*( C[11] + C41 );
422 C[10] += T*C[13];
423 C[ 6] = C30;
424 C[11] = C41;
425
426 C[15] += T*C[18];
427 C[16] += T54;
428 C[ 3] += T*(C[ 8] + C[15]);
429 C[ 4] += T*(C[12] + C[16]);
430 C[ 5] += T*(C[17] + C52);
431 C[17] = C52;
432 C[12] += T54;
433
434 C[ 7] += T*C[13];
435 C[ 8] += T*C[18];
436
437 C[21] += T*C[24];
438 C[22] += T*C[25];
439 C[23] += T*C[26];
440 C[28] += T*C[31];
441 C[29] += T*C[32];
442 C[30] += T*C[33];
443
444}
445
446void CbmKFSecondaryVertexFinder::AddTopoConstraint()
447{
448 if( !VParent ) return;
449
450 double m[3], *V;
451 {
452 m[0] = VParent->GetRefX();
453 m[1] = VParent->GetRefY();
454 m[2] = VParent->GetRefZ();
455 V = VParent->GetCovMatrix();
456 }
457
458 Extrapolate(-r0[7]);
459
460 {
461 double inf = 1.e8;
462 double xinf = r0[0]*inf;
463 double yinf = r0[1]*inf;
464 double zinf = r0[2]*inf;
465 C[0] += r0[0]*xinf;
466 C[1] += r0[0]*yinf;
467 C[2] += r0[1]*yinf;
468 C[3] += r0[0]*zinf;
469 C[4] += r0[1]*zinf;
470 C[5] += r0[2]*zinf;
471 C[28] = -2*xinf;
472 C[29] = -2*yinf;
473 C[30] = -2*zinf;
474 C[31] = 0;
475 C[32] = 0;
476 C[33] = 0;
477 C[34] = 0;
478 C[35] = inf;
479 }
480 r[7] = r0[7];
481
482 double Ai[6];
483 Ai[0] = C[4]*C[4] - C[2]*C[5];
484 Ai[1] = C[1]*C[5] - C[3]*C[4];
485 Ai[3] = C[2]*C[3] - C[1]*C[4];
486 double det = (C[0]*Ai[0] + C[1]*Ai[1] + C[3]*Ai[3]);
487 det = (det>1.e-8) ?1./det :0;
488 Ai[0] *= det;
489 Ai[1] *= det;
490 Ai[3] *= det;
491 Ai[2] = ( C[3]*C[3] - C[0]*C[5] )*det;
492 Ai[4] = ( C[0]*C[4] - C[1]*C[3] )*det;
493 Ai[5] = ( C[1]*C[1] - C[0]*C[2] )*det;
494
495 double B[5][3];
496 B[0][0] = C[ 6]*Ai[0] + C[ 7]*Ai[1] + C[ 8]*Ai[3];
497 B[0][1] = C[ 6]*Ai[1] + C[ 7]*Ai[2] + C[ 8]*Ai[4];
498 B[0][2] = C[ 6]*Ai[3] + C[ 7]*Ai[4] + C[ 8]*Ai[5];
499
500 B[1][0] = C[10]*Ai[0] + C[11]*Ai[1] + C[12]*Ai[3];
501 B[1][1] = C[10]*Ai[1] + C[11]*Ai[2] + C[12]*Ai[4];
502 B[1][2] = C[10]*Ai[3] + C[11]*Ai[4] + C[12]*Ai[5];
503
504 B[2][0] = C[15]*Ai[0] + C[16]*Ai[1] + C[17]*Ai[3];
505 B[2][1] = C[15]*Ai[1] + C[16]*Ai[2] + C[17]*Ai[4];
506 B[2][2] = C[15]*Ai[3] + C[16]*Ai[4] + C[17]*Ai[5];
507
508 B[3][0] = C[21]*Ai[0] + C[22]*Ai[1] + C[23]*Ai[3];
509 B[3][1] = C[21]*Ai[1] + C[22]*Ai[2] + C[23]*Ai[4];
510 B[3][2] = C[21]*Ai[3] + C[22]*Ai[4] + C[23]*Ai[5];
511
512 B[4][0] = C[28]*Ai[0] + C[29]*Ai[1] + C[30]*Ai[3];
513 B[4][1] = C[28]*Ai[1] + C[29]*Ai[2] + C[30]*Ai[4];
514 B[4][2] = C[28]*Ai[3] + C[29]*Ai[4] + C[30]*Ai[5];
515
516 double z[3] = {m[0]-r[0], m[1]-r[1], m[2]-r[2]};
517 r[0] = m[0];
518 r[1] = m[1];
519 r[2] = m[2];
520 r[3]+= B[0][0]*z[0] + B[0][1]*z[1] + B[0][2]*z[2];
521 r[4]+= B[1][0]*z[0] + B[1][1]*z[1] + B[1][2]*z[2];
522 r[5]+= B[2][0]*z[0] + B[2][1]*z[1] + B[2][2]*z[2];
523 r[6]+= B[3][0]*z[0] + B[3][1]*z[1] + B[3][2]*z[2];
524 r[7]+= B[4][0]*z[0] + B[4][1]*z[1] + B[4][2]*z[2];
525
526 {
527 double AV[6] = { C[0]-V[0], C[1]-V[1], C[2]-V[2],
528 C[3]-V[3], C[4]-V[4], C[5]-V[5] };
529 double AVi[6];
530 AVi[0] = AV[4]*AV[4] - AV[2]*AV[5];
531 AVi[1] = AV[1]*AV[5] - AV[3]*AV[4];
532 AVi[2] = AV[3]*AV[3] - AV[0]*AV[5] ;
533 AVi[3] = AV[2]*AV[3] - AV[1]*AV[4];
534 AVi[4] = AV[0]*AV[4] - AV[1]*AV[3] ;
535 AVi[5] = -123; // FIXME assing correct value
536 det = ( AV[0]*AVi[0] + AV[1]*AVi[1] + AV[3]*AVi[3] );
537 det = (det>1.e-8) ?1./det :0;
538 NDF += 2;
539 Chi2 += ( +(AVi[0]*z[0] + AVi[1]*z[1] + AVi[3]*z[2])*z[0]
540 +(AVi[1]*z[0] + AVi[2]*z[1] + AVi[4]*z[2])*z[1]
541 +(AVi[3]*z[0] + AVi[4]*z[1] + AVi[5]*z[2])*z[2] )*det;
542 }
543
544 double d0, d1, d2;
545 C[0] = V[0];
546 C[1] = V[1];
547 C[2] = V[2];
548 C[3] = V[3];
549 C[4] = V[4];
550 C[5] = V[5];
551
552 d0= B[0][0]*V[0] + B[0][1]*V[1] + B[0][2]*V[3] - C[ 6];
553 d1= B[0][0]*V[1] + B[0][1]*V[2] + B[0][2]*V[4] - C[ 7];
554 d2= B[0][0]*V[3] + B[0][1]*V[4] + B[0][2]*V[5] - C[ 8];
555
556 C[ 6]+= d0;
557 C[ 7]+= d1;
558 C[ 8]+= d2;
559 C[ 9]+= d0*B[0][0] + d1*B[0][1] + d2*B[0][2];
560
561 d0= B[1][0]*V[0] + B[1][1]*V[1] + B[1][2]*V[3] - C[10];
562 d1= B[1][0]*V[1] + B[1][1]*V[2] + B[1][2]*V[4] - C[11];
563 d2= B[1][0]*V[3] + B[1][1]*V[4] + B[1][2]*V[5] - C[12];
564
565 C[10]+= d0;
566 C[11]+= d1;
567 C[12]+= d2;
568 C[13]+= d0*B[0][0] + d1*B[0][1] + d2*B[0][2];
569 C[14]+= d0*B[1][0] + d1*B[1][1] + d2*B[1][2];
570
571 d0= B[2][0]*V[0] + B[2][1]*V[1] + B[2][2]*V[3] - C[15];
572 d1= B[2][0]*V[1] + B[2][1]*V[2] + B[2][2]*V[4] - C[16];
573 d2= B[2][0]*V[3] + B[2][1]*V[4] + B[2][2]*V[5] - C[17];
574
575 C[15]+= d0;
576 C[16]+= d1;
577 C[17]+= d2;
578 C[18]+= d0*B[0][0] + d1*B[0][1] + d2*B[0][2];
579 C[19]+= d0*B[1][0] + d1*B[1][1] + d2*B[1][2];
580 C[20]+= d0*B[2][0] + d1*B[2][1] + d2*B[2][2];
581
582 d0= B[3][0]*V[0] + B[3][1]*V[1] + B[3][2]*V[3] - C[21];
583 d1= B[3][0]*V[1] + B[3][1]*V[2] + B[3][2]*V[4] - C[22];
584 d2= B[3][0]*V[3] + B[3][1]*V[4] + B[3][2]*V[5] - C[23];
585
586 C[21]+= d0;
587 C[22]+= d1;
588 C[23]+= d2;
589 C[24]+= d0*B[0][0] + d1*B[0][1] + d2*B[0][2];
590 C[25]+= d0*B[1][0] + d1*B[1][1] + d2*B[1][2];
591 C[26]+= d0*B[2][0] + d1*B[2][1] + d2*B[2][2];
592 C[27]+= d0*B[3][0] + d1*B[3][1] + d2*B[3][2];
593
594 d0= B[4][0]*V[0] + B[4][1]*V[1] + B[4][2]*V[3] - C[28];
595 d1= B[4][0]*V[1] + B[4][1]*V[2] + B[4][2]*V[4] - C[29];
596 d2= B[4][0]*V[3] + B[4][1]*V[4] + B[4][2]*V[5] - C[30];
597
598 C[28]+= d0;
599 C[29]+= d1;
600 C[30]+= d2;
601 C[31]+= d0*B[0][0] + d1*B[0][1] + d2*B[0][2];
602 C[32]+= d0*B[1][0] + d1*B[1][1] + d2*B[1][2];
603 C[33]+= d0*B[2][0] + d1*B[2][1] + d2*B[2][2];
604 C[34]+= d0*B[3][0] + d1*B[3][1] + d2*B[3][2];
605 C[35]+= d0*B[4][0] + d1*B[4][1] + d2*B[4][2];
606
607 {
608 Extrapolate(r[7]);
609
610 double CH00 = C[0]+r0[0]*C[28];
611 double CH10 = C[1]+r0[0]*C[29];
612 double CH11 = C[2]+r0[1]*C[29];
613 double CH20 = C[3]+r0[0]*C[30];
614 double CH21 = C[4]+r0[1]*C[30];
615 double CH22 = C[5]+r0[2]*C[30];
616
617 C[28]+= r0[0]*C[35];
618 C[29]+= r0[1]*C[35];
619 C[30]+= r0[2]*C[35];
620
621 C[0] = CH00+r0[0]*C[28];
622 C[1] = CH10+r0[1]*C[28];
623 C[2] = CH11+r0[1]*C[29];
624 C[3] = CH20+r0[2]*C[28];
625 C[4] = CH21+r0[2]*C[29];
626 C[5] = CH22+r0[2]*C[30];
627
628 C[6]+= r0[0]*C[31];
629 C[7]+= r0[1]*C[31];
630 C[8]+= r0[2]*C[31];
631 C[10]+= r0[0]*C[32];
632 C[11]+= r0[1]*C[32];
633 C[12]+= r0[2]*C[32];
634 C[15]+= r0[0]*C[33];
635 C[16]+= r0[1]*C[33];
636 C[17]+= r0[2]*C[33];
637 C[21]+= r0[0]*C[34];
638 C[22]+= r0[1]*C[34];
639 C[23]+= r0[2]*C[34];
640
641 }
642}
friend F32vec4 sqrt(const F32vec4 &a)
Definition P4_F32vec4.h:34
int i
Definition P4_F32vec4.h:22
__m128 m
Definition P4_F32vec4.h:27
void SetTopoConstraint(CbmKFVertexInterface *Parent=0)
void SetTracks(std::vector< CbmKFTrackInterface * > &vTracks)
void GetMass(Double_t *M, Double_t *Error)
void SetApproximation(CbmKFVertexInterface *Guess=0)
void SetMassConstraint(Double_t MotherMass=-1)
void GetVertex(CbmKFVertexInterface &vtx)
void GetMotherTrack(Double_t T[], Double_t C[])
void AddTrack(CbmKFTrackInterface *Track)
Int_t Extrapolate(Double_t z, Double_t *QP0=0, Bool_t line=false)
Access to i-th hit.
Double_t GetMass()
Definition CbmKFTrack.h:55
Double_t * GetTrack()
Is it electron.
Definition CbmKFTrack.h:51
Double_t & GetRefChi2()
array[15] of covariance matrix
Definition CbmKFTrack.h:53
Int_t & GetRefNDF()
Chi^2 after fit.
Definition CbmKFTrack.h:54
Double_t * GetCovMatrix()
array[6] of track parameters(x,y,tx,ty,qp,z)
Definition CbmKFTrack.h:52
Double_t z
virtual Double_t & GetRefY()
virtual Double_t & GetRefX()
virtual Double_t & GetRefZ()
virtual Double_t & GetRefChi2()
Array[6] of covariance matrix.
virtual Int_t & GetRefNDF()
Chi^2 after fit.
void GetVertex(CbmVertex &v)
virtual Double_t * GetCovMatrix()
static CbmKF * Instance()
Definition CbmKF.h:35
std::vector< CbmKFTube > vTargets
Definition CbmKF.h:63