00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019 #include "KFParticleBase.h"
00020 #include "TMath.h"
00021 #include "Riostream.h"
00022 #include "TString.h"
00023 ClassImp(KFParticleBase);
00024
00025
00026 KFParticleBase::KFParticleBase() :fID(0), fQ(0), fNDF(-3), fChi2(0), fSFromDecay(0), fAtProductionVertex(0), fIsLinearized(0),
00027 fIdTruth(0), fQuality(0), fIdParentVx(0)
00028 {
00029
00030
00031 Initialize();
00032 }
00033
00034 void KFParticleBase::Initialize( const Double_t Param[], const Double_t Cov[], Int_t Charge, Double_t Mass )
00035 {
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049 for( Int_t i=0; i<6 ; i++ ) fP[i] = Param[i];
00050 for( Int_t i=0; i<21; i++ ) fC[i] = Cov[i];
00051
00052 Double_t energy = TMath::Sqrt( Mass*Mass + fP[3]*fP[3] + fP[4]*fP[4] + fP[5]*fP[5]);
00053 fP[6] = energy;
00054 fP[7] = 0;
00055 fQ = Charge;
00056 fNDF = 0;
00057 fChi2 = 0;
00058 fAtProductionVertex = 0;
00059 fIsLinearized = 0;
00060 fSFromDecay = 0;
00061
00062 Double_t energyInv = 1./energy;
00063 Double_t
00064 h0 = fP[3]*energyInv,
00065 h1 = fP[4]*energyInv,
00066 h2 = fP[5]*energyInv;
00067
00068 fC[21] = h0*fC[ 6] + h1*fC[10] + h2*fC[15];
00069 fC[22] = h0*fC[ 7] + h1*fC[11] + h2*fC[16];
00070 fC[23] = h0*fC[ 8] + h1*fC[12] + h2*fC[17];
00071 fC[24] = h0*fC[ 9] + h1*fC[13] + h2*fC[18];
00072 fC[25] = h0*fC[13] + h1*fC[14] + h2*fC[19];
00073 fC[26] = h0*fC[18] + h1*fC[19] + h2*fC[20];
00074 fC[27] = ( h0*h0*fC[ 9] + h1*h1*fC[14] + h2*h2*fC[20]
00075 + 2*(h0*h1*fC[13] + h0*h2*fC[18] + h1*h2*fC[19] ) );
00076 for( Int_t i=28; i<36; i++ ) fC[i] = 0;
00077 fC[35] = 1.;
00078 }
00079
00080 void KFParticleBase::Initialize()
00081 {
00082
00083
00084 for( Int_t i=0; i<8; i++) fP[i] = 0;
00085 for(Int_t i=0;i<36;++i) fC[i]=0.;
00086 fC[0] = fC[2] = fC[5] = 100.;
00087 fC[35] = 1.;
00088 fNDF = -3;
00089 fChi2 = 0.;
00090 fQ = 0;
00091 fSFromDecay = 0;
00092 fAtProductionVertex = 0;
00093 fVtxGuess[0]=fVtxGuess[1]=fVtxGuess[2]=0.;
00094 fIsLinearized = 0;
00095 }
00096
00097 void KFParticleBase::SetVtxGuess( Double_t x, Double_t y, Double_t z )
00098 {
00099
00100
00101 fVtxGuess[0] = x;
00102 fVtxGuess[1] = y;
00103 fVtxGuess[2] = z;
00104 fIsLinearized = 1;
00105 }
00106
00107
00108 Int_t KFParticleBase::GetMomentum( Double_t &p, Double_t &error ) const
00109 {
00110
00111
00112 Double_t x = fP[3];
00113 Double_t y = fP[4];
00114 Double_t z = fP[5];
00115 Double_t x2 = x*x;
00116 Double_t y2 = y*y;
00117 Double_t z2 = z*z;
00118 Double_t p2 = x2+y2+z2;
00119 p = TMath::Sqrt(p2);
00120 error = (x2*fC[9]+y2*fC[14]+z2*fC[20] + 2*(x*y*fC[13]+x*z*fC[18]+y*z*fC[19]) );
00121 if( error>0 && p>1.e-4 ){
00122 error = TMath::Sqrt(error)/p;
00123 return 0;
00124 }
00125 return 1;
00126 }
00127
00128 Int_t KFParticleBase::GetPt( Double_t &pt, Double_t &error ) const
00129 {
00130
00131
00132 Double_t px = fP[3];
00133 Double_t py = fP[4];
00134 Double_t px2 = px*px;
00135 Double_t py2 = py*py;
00136 Double_t pt2 = px2+py2;
00137 pt = TMath::Sqrt(pt2);
00138 error = (px2*fC[9] + py2*fC[14] + 2*px*py*fC[13] );
00139 if( error>0 && pt>1.e-4 ){
00140 error = TMath::Sqrt(error)/pt;
00141 return 0;
00142 }
00143 error = 1.e10;
00144 return 1;
00145 }
00146
00147 Int_t KFParticleBase::GetEta( Double_t &eta, Double_t &error ) const
00148 {
00149
00150
00151 Double_t px = fP[3];
00152 Double_t py = fP[4];
00153 Double_t pz = fP[5];
00154 Double_t pt2 = px*px + py*py;
00155 Double_t p2 = pt2 + pz*pz;
00156 Double_t p = TMath::Sqrt(p2);
00157 Double_t a = p + pz;
00158 Double_t b = p - pz;
00159 eta = 1.e10;
00160 if( b > 1.e-8 ){
00161 Double_t c = a/b;
00162 if( c>1.e-8 ) eta = 0.5*TMath::Log(c);
00163 }
00164 Double_t h3 = -px*pz;
00165 Double_t h4 = -py*pz;
00166 Double_t pt4 = pt2*pt2;
00167 Double_t p2pt4 = p2*pt4;
00168 error = (h3*h3*fC[9] + h4*h4*fC[14] + pt4*fC[20] + 2*( h3*(h4*fC[13] + fC[18]*pt2) + pt2*h4*fC[19] ) );
00169
00170 if( error>0 && p2pt4>1.e-10 ){
00171 error = TMath::Sqrt(error/p2pt4);
00172 return 0;
00173 }
00174
00175 error = 1.e10;
00176 return 1;
00177 }
00178
00179 Int_t KFParticleBase::GetPhi( Double_t &phi, Double_t &error ) const
00180 {
00181
00182
00183 Double_t px = fP[3];
00184 Double_t py = fP[4];
00185 Double_t px2 = px*px;
00186 Double_t py2 = py*py;
00187 Double_t pt2 = px2 + py2;
00188 phi = TMath::ATan2(py,px);
00189 error = (py2*fC[9] + px2*fC[14] - 2*px*py*fC[13] );
00190 if( error>0 && pt2>1.e-4 ){
00191 error = TMath::Sqrt(error)/pt2;
00192 return 0;
00193 }
00194 error = 1.e10;
00195 return 1;
00196 }
00197
00198 Int_t KFParticleBase::GetR( Double_t &r, Double_t &error ) const
00199 {
00200
00201
00202 Double_t x = fP[0];
00203 Double_t y = fP[1];
00204 Double_t x2 = x*x;
00205 Double_t y2 = y*y;
00206 r = TMath::Sqrt(x2 + y2);
00207 error = (x2*fC[0] + y2*fC[2] - 2*x*y*fC[1] );
00208 if( error>0 && r>1.e-4 ){
00209 error = TMath::Sqrt(error)/r;
00210 return 0;
00211 }
00212 error = 1.e10;
00213 return 1;
00214 }
00215
00216 Int_t KFParticleBase::GetMass( Double_t &m, Double_t &error ) const
00217 {
00218
00219
00220
00221
00222 Double_t s = ( fP[3]*fP[3]*fC[9] + fP[4]*fP[4]*fC[14] + fP[5]*fP[5]*fC[20]
00223 + fP[6]*fP[6]*fC[27]
00224 +2*( + fP[3]*fP[4]*fC[13] + fP[5]*(fP[3]*fC[18] + fP[4]*fC[19])
00225 - fP[6]*( fP[3]*fC[24] + fP[4]*fC[25] + fP[5]*fC[26] ) )
00226 );
00227 Double_t m2 = TMath::Abs(fP[6]*fP[6] - fP[3]*fP[3] - fP[4]*fP[4] - fP[5]*fP[5]);
00228 m = TMath::Sqrt(m2);
00229 if( m>1.e-10 ){
00230 if( s>=0 ){
00231 error = TMath::Sqrt(s)/m;
00232 return 0;
00233 }
00234 }
00235 error = 1.e20;
00236 return 1;
00237 }
00238
00239
00240 Int_t KFParticleBase::GetDecayLength( Double_t &l, Double_t &error ) const
00241 {
00242
00243
00244 Double_t x = fP[3];
00245 Double_t y = fP[4];
00246 Double_t z = fP[5];
00247 Double_t t = fP[7];
00248 Double_t x2 = x*x;
00249 Double_t y2 = y*y;
00250 Double_t z2 = z*z;
00251 Double_t p2 = x2+y2+z2;
00252 l = t*TMath::Sqrt(p2);
00253 if( p2>1.e-4){
00254 error = p2*fC[35] + t*t/p2*(x2*fC[9]+y2*fC[14]+z2*fC[20]
00255 + 2*(x*y*fC[13]+x*z*fC[18]+y*z*fC[19]) )
00256 + 2*t*(x*fC[31]+y*fC[32]+z*fC[33]);
00257 error = TMath::Sqrt(TMath::Abs(error));
00258 return 0;
00259 }
00260 error = 1.e20;
00261 return 1;
00262 }
00263
00264 Int_t KFParticleBase::GetDecayLengthXY( Double_t &l, Double_t &error ) const
00265 {
00266
00267
00268 Double_t x = fP[3];
00269 Double_t y = fP[4];
00270 Double_t t = fP[7];
00271 Double_t x2 = x*x;
00272 Double_t y2 = y*y;
00273 Double_t pt2 = x2+y2;
00274 l = t*TMath::Sqrt(pt2);
00275 if( pt2>1.e-4){
00276 error = pt2*fC[35] + t*t/pt2*(x2*fC[9]+y2*fC[14] + 2*x*y*fC[13] )
00277 + 2*t*(x*fC[31]+y*fC[32]);
00278 error = TMath::Sqrt(TMath::Abs(error));
00279 return 0;
00280 }
00281 error = 1.e20;
00282 return 1;
00283 }
00284
00285
00286 Int_t KFParticleBase::GetLifeTime( Double_t &tauC, Double_t &error ) const
00287 {
00288
00289
00290 Double_t m, dm;
00291 GetMass( m, dm );
00292 Double_t cTM = (-fP[3]*fC[31] - fP[4]*fC[32] - fP[5]*fC[33] + fP[6]*fC[34]);
00293 tauC = fP[7]*m;
00294 error = m*m*fC[35] + 2*fP[7]*cTM + fP[7]*fP[7]*dm*dm;
00295 if( error > 0 ){
00296 error = TMath::Sqrt( error );
00297 return 0;
00298 }
00299 error = 1.e20;
00300 return 1;
00301 }
00302
00303
00304 void KFParticleBase::operator +=( const KFParticleBase &Daughter )
00305 {
00306
00307
00308 AddDaughter( Daughter );
00309 }
00310
00311 Double_t KFParticleBase::GetSCorrection( const Double_t Part[], const Double_t XYZ[] )
00312 {
00313
00314
00315 Double_t d[3] = { XYZ[0]-Part[0], XYZ[1]-Part[1], XYZ[2]-Part[2] };
00316 Double_t p2 = Part[3]*Part[3]+Part[4]*Part[4]+Part[5]*Part[5];
00317 Double_t sigmaS = (p2>1.e-4) ? ( .1+3.*TMath::Sqrt( d[0]*d[0]+d[1]*d[1]+d[2]*d[2]) )/TMath::Sqrt(p2) : 1.;
00318 return sigmaS;
00319 }
00320
00321 void KFParticleBase::GetMeasurement( const Double_t XYZ[], Double_t m[], Double_t V[] ) const
00322 {
00323
00324
00325 Double_t b[3];
00326 GetFieldValue( XYZ, b );
00327 const Double_t kCLight = 0.000299792458;
00328 b[0]*=kCLight; b[1]*=kCLight; b[2]*=kCLight;
00329
00330 Transport( GetDStoPoint(XYZ), m, V );
00331
00332 Double_t sigmaS = GetSCorrection( m, XYZ );
00333
00334 Double_t h[6];
00335
00336 h[0] = m[3]*sigmaS;
00337 h[1] = m[4]*sigmaS;
00338 h[2] = m[5]*sigmaS;
00339 h[3] = ( h[1]*b[2]-h[2]*b[1] )*GetQ();
00340 h[4] = ( h[2]*b[0]-h[0]*b[2] )*GetQ();
00341 h[5] = ( h[0]*b[1]-h[1]*b[0] )*GetQ();
00342
00343 V[ 0]+= h[0]*h[0];
00344 V[ 1]+= h[1]*h[0];
00345 V[ 2]+= h[1]*h[1];
00346 V[ 3]+= h[2]*h[0];
00347 V[ 4]+= h[2]*h[1];
00348 V[ 5]+= h[2]*h[2];
00349
00350 V[ 6]+= h[3]*h[0];
00351 V[ 7]+= h[3]*h[1];
00352 V[ 8]+= h[3]*h[2];
00353 V[ 9]+= h[3]*h[3];
00354
00355 V[10]+= h[4]*h[0];
00356 V[11]+= h[4]*h[1];
00357 V[12]+= h[4]*h[2];
00358 V[13]+= h[4]*h[3];
00359 V[14]+= h[4]*h[4];
00360
00361 V[15]+= h[5]*h[0];
00362 V[16]+= h[5]*h[1];
00363 V[17]+= h[5]*h[2];
00364 V[18]+= h[5]*h[3];
00365 V[19]+= h[5]*h[4];
00366 V[20]+= h[5]*h[5];
00367 }
00368
00369
00370 void KFParticleBase::AddDaughter( const KFParticleBase &Daughter )
00371 {
00372
00373
00374 if( fNDF<-1 ){
00375 fNDF = -1;
00376 fQ = Daughter.GetQ();
00377 for( Int_t i=0; i<7; i++) fP[i] = Daughter.fP[i];
00378 for( Int_t i=0; i<28; i++) fC[i] = Daughter.fC[i];
00379 fSFromDecay = 0;
00380 return;
00381 }
00382
00383 TransportToDecayVertex();
00384
00385 Double_t b[3];
00386 Int_t maxIter = 1;
00387
00388 if( !fIsLinearized ){
00389 if( fNDF==-1 ){
00390 Double_t ds, ds1;
00391 GetDStoParticle(Daughter, ds, ds1);
00392 TransportToDS( ds );
00393 Double_t m[8];
00394 Double_t mCd[36];
00395 Daughter.Transport( ds1, m, mCd );
00396 fVtxGuess[0] = .5*( fP[0] + m[0] );
00397 fVtxGuess[1] = .5*( fP[1] + m[1] );
00398 fVtxGuess[2] = .5*( fP[2] + m[2] );
00399 } else {
00400 fVtxGuess[0] = fP[0];
00401 fVtxGuess[1] = fP[1];
00402 fVtxGuess[2] = fP[2];
00403 }
00404 maxIter = 3;
00405 }
00406
00407 for( Int_t iter=0; iter<maxIter; iter++ ){
00408
00409 {
00410 GetFieldValue( fVtxGuess, b );
00411 const Double_t kCLight = 0.000299792458;
00412 b[0]*=kCLight; b[1]*=kCLight; b[2]*=kCLight;
00413 }
00414
00415 Double_t *ffP = fP, *ffC = fC, tmpP[8], tmpC[36];
00416 if( fNDF==-1 ){
00417 GetMeasurement( fVtxGuess, tmpP, tmpC );
00418 ffP = tmpP;
00419 ffC = tmpC;
00420 }
00421
00422 Double_t m[8], mV[36];
00423
00424 if( Daughter.fC[35]>0 ){
00425 Daughter.GetMeasurement( fVtxGuess, m, mV );
00426 } else {
00427 for( Int_t i=0; i<8; i++ ) m[i] = Daughter.fP[i];
00428 for( Int_t i=0; i<36; i++ ) mV[i] = Daughter.fC[i];
00429 }
00430
00431
00432
00433 Double_t mS[6];
00434 {
00435 Double_t mSi[6] = { ffC[0]+mV[0],
00436 ffC[1]+mV[1], ffC[2]+mV[2],
00437 ffC[3]+mV[3], ffC[4]+mV[4], ffC[5]+mV[5] };
00438
00439 mS[0] = mSi[2]*mSi[5] - mSi[4]*mSi[4];
00440 mS[1] = mSi[3]*mSi[4] - mSi[1]*mSi[5];
00441 mS[2] = mSi[0]*mSi[5] - mSi[3]*mSi[3];
00442 mS[3] = mSi[1]*mSi[4] - mSi[2]*mSi[3];
00443 mS[4] = mSi[1]*mSi[3] - mSi[0]*mSi[4];
00444 mS[5] = mSi[0]*mSi[2] - mSi[1]*mSi[1];
00445
00446 Double_t s = ( mSi[0]*mS[0] + mSi[1]*mS[1] + mSi[3]*mS[3] );
00447
00448 s = ( s > 1.E-20 ) ?1./s :0;
00449 mS[0]*=s;
00450 mS[1]*=s;
00451 mS[2]*=s;
00452 mS[3]*=s;
00453 mS[4]*=s;
00454 mS[5]*=s;
00455 }
00456
00457
00458
00459 Double_t zeta[3] = { m[0]-ffP[0], m[1]-ffP[1], m[2]-ffP[2] };
00460
00461
00462
00463
00464 Double_t mCHt0[7], mCHt1[7], mCHt2[7];
00465
00466 mCHt0[0]=ffC[ 0] ; mCHt1[0]=ffC[ 1] ; mCHt2[0]=ffC[ 3] ;
00467 mCHt0[1]=ffC[ 1] ; mCHt1[1]=ffC[ 2] ; mCHt2[1]=ffC[ 4] ;
00468 mCHt0[2]=ffC[ 3] ; mCHt1[2]=ffC[ 4] ; mCHt2[2]=ffC[ 5] ;
00469 mCHt0[3]=ffC[ 6]-mV[ 6]; mCHt1[3]=ffC[ 7]-mV[ 7]; mCHt2[3]=ffC[ 8]-mV[ 8];
00470 mCHt0[4]=ffC[10]-mV[10]; mCHt1[4]=ffC[11]-mV[11]; mCHt2[4]=ffC[12]-mV[12];
00471 mCHt0[5]=ffC[15]-mV[15]; mCHt1[5]=ffC[16]-mV[16]; mCHt2[5]=ffC[17]-mV[17];
00472 mCHt0[6]=ffC[21]-mV[21]; mCHt1[6]=ffC[22]-mV[22]; mCHt2[6]=ffC[23]-mV[23];
00473
00474
00475
00476 Double_t k0[7], k1[7], k2[7];
00477
00478 for(Int_t i=0;i<7;++i){
00479 k0[i] = mCHt0[i]*mS[0] + mCHt1[i]*mS[1] + mCHt2[i]*mS[3];
00480 k1[i] = mCHt0[i]*mS[1] + mCHt1[i]*mS[2] + mCHt2[i]*mS[4];
00481 k2[i] = mCHt0[i]*mS[3] + mCHt1[i]*mS[4] + mCHt2[i]*mS[5];
00482 }
00483
00484
00485
00486 if( iter<maxIter-1 ){
00487 for(Int_t i=0; i<3; ++i)
00488 fVtxGuess[i]= ffP[i] + k0[i]*zeta[0]+k1[i]*zeta[1]+k2[i]*zeta[2];
00489 continue;
00490 }
00491
00492
00493
00494
00495
00496 ffP[ 3] += m[ 3];
00497 ffP[ 4] += m[ 4];
00498 ffP[ 5] += m[ 5];
00499 ffP[ 6] += m[ 6];
00500
00501 ffC[ 9] += mV[ 9];
00502 ffC[13] += mV[13];
00503 ffC[14] += mV[14];
00504 ffC[18] += mV[18];
00505 ffC[19] += mV[19];
00506 ffC[20] += mV[20];
00507 ffC[24] += mV[24];
00508 ffC[25] += mV[25];
00509 ffC[26] += mV[26];
00510 ffC[27] += mV[27];
00511
00512
00513
00514
00515 for(Int_t i=0;i<7;++i)
00516 fP[i] = ffP[i] + k0[i]*zeta[0] + k1[i]*zeta[1] + k2[i]*zeta[2];
00517
00518
00519
00520 for(Int_t i=0, k=0;i<7;++i){
00521 for(Int_t j=0;j<=i;++j,++k) {
00522 fC[k] = ffC[k] - (k0[i]*mCHt0[j] + k1[i]*mCHt1[j] + k2[i]*mCHt2[j] );
00523 }
00524 }
00525
00526
00527
00528 fNDF += 2;
00529 fQ += Daughter.GetQ();
00530 fSFromDecay = 0;
00531 fChi2 += (mS[0]*zeta[0] + mS[1]*zeta[1] + mS[3]*zeta[2])*zeta[0]
00532 + (mS[1]*zeta[0] + mS[2]*zeta[1] + mS[4]*zeta[2])*zeta[1]
00533 + (mS[3]*zeta[0] + mS[4]*zeta[1] + mS[5]*zeta[2])*zeta[2];
00534
00535 }
00536 }
00537
00538
00539 void KFParticleBase::SetProductionVertex( const KFParticleBase &Vtx )
00540 {
00541
00542
00543 const Double_t *m = Vtx.fP, *mV = Vtx.fC;
00544
00545 Bool_t noS = ( fC[35]<=0 );
00546
00547 if( noS ){
00548 TransportToDecayVertex();
00549 fP[7] = 0;
00550 fC[28] = 0; fC[29] = 0; fC[30] = 0; fC[31] = 0; fC[32] = 0; fC[33] = 0; fC[35] = 0; fC[35] = 0;
00551 } else {
00552 TransportToDS( GetDStoPoint( m ) );
00553 fP[7] = -fSFromDecay;
00554 Convert(1);
00555 }
00556
00557 Double_t mAi[6];
00558
00559 InvertSym3( fC, mAi );
00560
00561 Double_t mB[5][3];
00562
00563 mB[0][0] = fC[ 6]*mAi[0] + fC[ 7]*mAi[1] + fC[ 8]*mAi[3];
00564 mB[0][1] = fC[ 6]*mAi[1] + fC[ 7]*mAi[2] + fC[ 8]*mAi[4];
00565 mB[0][2] = fC[ 6]*mAi[3] + fC[ 7]*mAi[4] + fC[ 8]*mAi[5];
00566
00567 mB[1][0] = fC[10]*mAi[0] + fC[11]*mAi[1] + fC[12]*mAi[3];
00568 mB[1][1] = fC[10]*mAi[1] + fC[11]*mAi[2] + fC[12]*mAi[4];
00569 mB[1][2] = fC[10]*mAi[3] + fC[11]*mAi[4] + fC[12]*mAi[5];
00570
00571 mB[2][0] = fC[15]*mAi[0] + fC[16]*mAi[1] + fC[17]*mAi[3];
00572 mB[2][1] = fC[15]*mAi[1] + fC[16]*mAi[2] + fC[17]*mAi[4];
00573 mB[2][2] = fC[15]*mAi[3] + fC[16]*mAi[4] + fC[17]*mAi[5];
00574
00575 mB[3][0] = fC[21]*mAi[0] + fC[22]*mAi[1] + fC[23]*mAi[3];
00576 mB[3][1] = fC[21]*mAi[1] + fC[22]*mAi[2] + fC[23]*mAi[4];
00577 mB[3][2] = fC[21]*mAi[3] + fC[22]*mAi[4] + fC[23]*mAi[5];
00578
00579 mB[4][0] = fC[28]*mAi[0] + fC[29]*mAi[1] + fC[30]*mAi[3];
00580 mB[4][1] = fC[28]*mAi[1] + fC[29]*mAi[2] + fC[30]*mAi[4];
00581 mB[4][2] = fC[28]*mAi[3] + fC[29]*mAi[4] + fC[30]*mAi[5];
00582
00583 Double_t z[3] = { m[0]-fP[0], m[1]-fP[1], m[2]-fP[2] };
00584
00585 {
00586 Double_t mAVi[6] = { fC[0]-mV[0], fC[1]-mV[1], fC[2]-mV[2],
00587 fC[3]-mV[3], fC[4]-mV[4], fC[5]-mV[5] };
00588
00589 if( !InvertSym3( mAVi, mAVi ) ){
00590
00591 Double_t dChi2 = ( +(mAVi[0]*z[0] + mAVi[1]*z[1] + mAVi[3]*z[2])*z[0]
00592 +(mAVi[1]*z[0] + mAVi[2]*z[1] + mAVi[4]*z[2])*z[1]
00593 +(mAVi[3]*z[0] + mAVi[4]*z[1] + mAVi[5]*z[2])*z[2] );
00594
00595
00596
00597
00598 fChi2+= TMath::Abs( dChi2 );
00599 }
00600 fNDF += 2;
00601 }
00602
00603 fP[0] = m[0];
00604 fP[1] = m[1];
00605 fP[2] = m[2];
00606 fP[3]+= mB[0][0]*z[0] + mB[0][1]*z[1] + mB[0][2]*z[2];
00607 fP[4]+= mB[1][0]*z[0] + mB[1][1]*z[1] + mB[1][2]*z[2];
00608 fP[5]+= mB[2][0]*z[0] + mB[2][1]*z[1] + mB[2][2]*z[2];
00609 fP[6]+= mB[3][0]*z[0] + mB[3][1]*z[1] + mB[3][2]*z[2];
00610 fP[7]+= mB[4][0]*z[0] + mB[4][1]*z[1] + mB[4][2]*z[2];
00611
00612 Double_t d0, d1, d2;
00613
00614 fC[0] = mV[0];
00615 fC[1] = mV[1];
00616 fC[2] = mV[2];
00617 fC[3] = mV[3];
00618 fC[4] = mV[4];
00619 fC[5] = mV[5];
00620
00621 d0= mB[0][0]*mV[0] + mB[0][1]*mV[1] + mB[0][2]*mV[3] - fC[ 6];
00622 d1= mB[0][0]*mV[1] + mB[0][1]*mV[2] + mB[0][2]*mV[4] - fC[ 7];
00623 d2= mB[0][0]*mV[3] + mB[0][1]*mV[4] + mB[0][2]*mV[5] - fC[ 8];
00624
00625 fC[ 6]+= d0;
00626 fC[ 7]+= d1;
00627 fC[ 8]+= d2;
00628 fC[ 9]+= d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
00629
00630 d0= mB[1][0]*mV[0] + mB[1][1]*mV[1] + mB[1][2]*mV[3] - fC[10];
00631 d1= mB[1][0]*mV[1] + mB[1][1]*mV[2] + mB[1][2]*mV[4] - fC[11];
00632 d2= mB[1][0]*mV[3] + mB[1][1]*mV[4] + mB[1][2]*mV[5] - fC[12];
00633
00634 fC[10]+= d0;
00635 fC[11]+= d1;
00636 fC[12]+= d2;
00637 fC[13]+= d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
00638 fC[14]+= d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
00639
00640 d0= mB[2][0]*mV[0] + mB[2][1]*mV[1] + mB[2][2]*mV[3] - fC[15];
00641 d1= mB[2][0]*mV[1] + mB[2][1]*mV[2] + mB[2][2]*mV[4] - fC[16];
00642 d2= mB[2][0]*mV[3] + mB[2][1]*mV[4] + mB[2][2]*mV[5] - fC[17];
00643
00644 fC[15]+= d0;
00645 fC[16]+= d1;
00646 fC[17]+= d2;
00647 fC[18]+= d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
00648 fC[19]+= d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
00649 fC[20]+= d0*mB[2][0] + d1*mB[2][1] + d2*mB[2][2];
00650
00651 d0= mB[3][0]*mV[0] + mB[3][1]*mV[1] + mB[3][2]*mV[3] - fC[21];
00652 d1= mB[3][0]*mV[1] + mB[3][1]*mV[2] + mB[3][2]*mV[4] - fC[22];
00653 d2= mB[3][0]*mV[3] + mB[3][1]*mV[4] + mB[3][2]*mV[5] - fC[23];
00654
00655 fC[21]+= d0;
00656 fC[22]+= d1;
00657 fC[23]+= d2;
00658 fC[24]+= d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
00659 fC[25]+= d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
00660 fC[26]+= d0*mB[2][0] + d1*mB[2][1] + d2*mB[2][2];
00661 fC[27]+= d0*mB[3][0] + d1*mB[3][1] + d2*mB[3][2];
00662
00663 d0= mB[4][0]*mV[0] + mB[4][1]*mV[1] + mB[4][2]*mV[3] - fC[28];
00664 d1= mB[4][0]*mV[1] + mB[4][1]*mV[2] + mB[4][2]*mV[4] - fC[29];
00665 d2= mB[4][0]*mV[3] + mB[4][1]*mV[4] + mB[4][2]*mV[5] - fC[30];
00666
00667 fC[28]+= d0;
00668 fC[29]+= d1;
00669 fC[30]+= d2;
00670 fC[31]+= d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
00671 fC[32]+= d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
00672 fC[33]+= d0*mB[2][0] + d1*mB[2][1] + d2*mB[2][2];
00673 fC[34]+= d0*mB[3][0] + d1*mB[3][1] + d2*mB[3][2];
00674 fC[35]+= d0*mB[4][0] + d1*mB[4][1] + d2*mB[4][2];
00675
00676 if( noS ){
00677 fP[7] = 0;
00678 fC[28] = 0;fC[29] = 0;fC[30] = 0;fC[31] = 0;fC[32] = 0;fC[33] = 0;fC[35] = 0;fC[35] = 0;
00679 } else {
00680 TransportToDS( fP[7] );
00681 Convert(0);
00682 }
00683
00684 fSFromDecay = 0;
00685 }
00686
00687
00688
00689 void KFParticleBase::SetMassConstraint( Double_t Mass, Double_t SigmaMass )
00690 {
00691
00692
00693 Double_t m2 = Mass*Mass;
00694 Double_t s2 = m2*SigmaMass*SigmaMass;
00695
00696 Double_t p2 = fP[3]*fP[3] + fP[4]*fP[4] + fP[5]*fP[5];
00697 Double_t e0 = TMath::Sqrt(m2+p2);
00698
00699 Double_t mH[8];
00700 mH[0] = mH[1] = mH[2] = 0.;
00701 mH[3] = -2*fP[3];
00702 mH[4] = -2*fP[4];
00703 mH[5] = -2*fP[5];
00704 mH[6] = 2*fP[6];
00705 mH[7] = 0;
00706
00707 Double_t zeta = e0*e0 - e0*fP[6];
00708 zeta = m2 - (fP[6]*fP[6]-p2);
00709
00710 Double_t mCHt[8], s2_est=0;
00711 for( Int_t i=0; i<8; ++i ){
00712 mCHt[i] = 0.0;
00713 for (Int_t j=0;j<8;++j) mCHt[i] += Cij(i,j)*mH[j];
00714 s2_est += mH[i]*mCHt[i];
00715 }
00716
00717 if( s2_est<1.e-20 ) return;
00718
00719
00720 Double_t w2 = 1./( s2 + s2_est );
00721 fChi2 += zeta*zeta*w2;
00722 fNDF += 1;
00723 for( Int_t i=0, ii=0; i<8; ++i ){
00724 Double_t ki = mCHt[i]*w2;
00725 fP[i]+= ki*zeta;
00726 for(Int_t j=0;j<=i;++j) fC[ii++] -= ki*mCHt[j];
00727 }
00728 }
00729
00730
00731 void KFParticleBase::SetNoDecayLength()
00732 {
00733
00734
00735 TransportToDecayVertex();
00736
00737 Double_t h[8];
00738 h[0] = h[1] = h[2] = h[3] = h[4] = h[5] = h[6] = 0;
00739 h[7] = 1;
00740
00741 Double_t zeta = 0 - fP[7];
00742 for(Int_t i=0;i<8;++i) zeta -= h[i]*(fP[i]-fP[i]);
00743
00744 Double_t s = fC[35];
00745 if( s>1.e-20 ){
00746 s = 1./s;
00747 fChi2 += zeta*zeta*s;
00748 fNDF += 1;
00749 for( Int_t i=0, ii=0; i<7; ++i ){
00750 Double_t ki = fC[28+i]*s;
00751 fP[i]+= ki*zeta;
00752 for(Int_t j=0;j<=i;++j) fC[ii++] -= ki*fC[28+j];
00753 }
00754 }
00755 fP[7] = 0;
00756 fC[28] = 0;fC[29] = 0;fC[30] = 0;fC[31] = 0;fC[32] = 0;fC[33] = 0;fC[35] = 0;fC[35] = 0;
00757 }
00758
00759
00760 void KFParticleBase::Construct( const KFParticleBase* vDaughters[], Int_t NDaughters,
00761 const KFParticleBase *Parent, Double_t Mass, Bool_t IsConstrained )
00762 {
00763
00764
00765 Int_t maxIter = 1;
00766 Bool_t wasLinearized = fIsLinearized;
00767 if( !fIsLinearized || IsConstrained ){
00768
00769 fVtxGuess[0] = GetX();
00770 fVtxGuess[1] = GetY();
00771 fVtxGuess[2] = GetZ();
00772 fIsLinearized = 1;
00773 maxIter = 3;
00774 }
00775
00776 Double_t constraintC[6];
00777
00778 if( IsConstrained ){
00779 for(Int_t i=0;i<6;++i) constraintC[i]=fC[i];
00780 } else {
00781 for(Int_t i=0;i<6;++i) constraintC[i]=0.;
00782 constraintC[0] = constraintC[2] = constraintC[5] = 100.;
00783 }
00784
00785
00786 for( Int_t iter=0; iter<maxIter; iter++ ){
00787 fAtProductionVertex = 0;
00788 fSFromDecay = 0;
00789 fP[0] = fVtxGuess[0];
00790 fP[1] = fVtxGuess[1];
00791 fP[2] = fVtxGuess[2];
00792 fP[3] = 0;
00793 fP[4] = 0;
00794 fP[5] = 0;
00795 fP[6] = 0;
00796 fP[7] = 0;
00797
00798 for(Int_t i=0;i<6; ++i) fC[i]=constraintC[i];
00799 for(Int_t i=6;i<36;++i) fC[i]=0.;
00800 fC[35] = 1.;
00801
00802 fNDF = IsConstrained ?0 :-3;
00803 fChi2 = 0.;
00804 fQ = 0;
00805
00806 for( Int_t itr =0; itr<NDaughters; itr++ ){
00807 AddDaughter( *vDaughters[itr] );
00808 }
00809 if( iter<maxIter-1){
00810 for( Int_t i=0; i<3; i++ ) fVtxGuess[i] = fP[i];
00811 }
00812 }
00813 fIsLinearized = wasLinearized;
00814
00815 if( Mass>=0 ) SetMassConstraint( Mass );
00816 if( Parent ) SetProductionVertex( *Parent );
00817 }
00818
00819
00820 void KFParticleBase::Convert( Bool_t ToProduction )
00821 {
00822
00823
00824
00825
00826 Double_t fld[3];
00827 {
00828 GetFieldValue( fP, fld );
00829 const Double_t kCLight = fQ*0.000299792458;
00830 fld[0]*=kCLight; fld[1]*=kCLight; fld[2]*=kCLight;
00831 }
00832
00833 Double_t h[6];
00834
00835 h[0] = fP[3];
00836 h[1] = fP[4];
00837 h[2] = fP[5];
00838 if( ToProduction ){ h[0]=-h[0]; h[1]=-h[1]; h[2]=-h[2]; }
00839 h[3] = h[1]*fld[2]-h[2]*fld[1];
00840 h[4] = h[2]*fld[0]-h[0]*fld[2];
00841 h[5] = h[0]*fld[1]-h[1]*fld[0];
00842
00843 Double_t c;
00844
00845 c = fC[28]+h[0]*fC[35];
00846 fC[ 0]+= h[0]*(c+fC[28]);
00847 fC[28] = c;
00848
00849 fC[ 1]+= h[1]*fC[28] + h[0]*fC[29];
00850 c = fC[29]+h[1]*fC[35];
00851 fC[ 2]+= h[1]*(c+fC[29]);
00852 fC[29] = c;
00853
00854 fC[ 3]+= h[2]*fC[28] + h[0]*fC[30];
00855 fC[ 4]+= h[2]*fC[29] + h[1]*fC[30];
00856 c = fC[30]+h[2]*fC[35];
00857 fC[ 5]+= h[2]*(c+fC[30]);
00858 fC[30] = c;
00859
00860 fC[ 6]+= h[3]*fC[28] + h[0]*fC[31];
00861 fC[ 7]+= h[3]*fC[29] + h[1]*fC[31];
00862 fC[ 8]+= h[3]*fC[30] + h[2]*fC[31];
00863 c = fC[31]+h[3]*fC[35];
00864 fC[ 9]+= h[3]*(c+fC[31]);
00865 fC[31] = c;
00866
00867 fC[10]+= h[4]*fC[28] + h[0]*fC[32];
00868 fC[11]+= h[4]*fC[29] + h[1]*fC[32];
00869 fC[12]+= h[4]*fC[30] + h[2]*fC[32];
00870 fC[13]+= h[4]*fC[31] + h[3]*fC[32];
00871 c = fC[32]+h[4]*fC[35];
00872 fC[14]+= h[4]*(c+fC[32]);
00873 fC[32] = c;
00874
00875 fC[15]+= h[5]*fC[28] + h[0]*fC[33];
00876 fC[16]+= h[5]*fC[29] + h[1]*fC[33];
00877 fC[17]+= h[5]*fC[30] + h[2]*fC[33];
00878 fC[18]+= h[5]*fC[31] + h[3]*fC[33];
00879 fC[19]+= h[5]*fC[32] + h[4]*fC[33];
00880 c = fC[33]+h[5]*fC[35];
00881 fC[20]+= h[5]*(c+fC[33]);
00882 fC[33] = c;
00883
00884 fC[21]+= h[0]*fC[34];
00885 fC[22]+= h[1]*fC[34];
00886 fC[23]+= h[2]*fC[34];
00887 fC[24]+= h[3]*fC[34];
00888 fC[25]+= h[4]*fC[34];
00889 fC[26]+= h[5]*fC[34];
00890 }
00891
00892
00893 void KFParticleBase::TransportToDecayVertex()
00894 {
00895
00896
00897 if( fSFromDecay != 0 ) TransportToDS( -fSFromDecay );
00898 if( fAtProductionVertex ) Convert(0);
00899 fAtProductionVertex = 0;
00900 }
00901
00902 void KFParticleBase::TransportToProductionVertex()
00903 {
00904
00905
00906 if( fSFromDecay != -fP[7] ) TransportToDS( -fSFromDecay-fP[7] );
00907 if( !fAtProductionVertex ) Convert( 1 );
00908 fAtProductionVertex = 1;
00909 }
00910
00911
00912 void KFParticleBase::TransportToDS( Double_t dS )
00913 {
00914
00915
00916 Transport( dS, fP, fC );
00917 fSFromDecay+= dS;
00918 }
00919
00920
00921 Double_t KFParticleBase::GetDStoPointLine( const Double_t xyz[] ) const
00922 {
00923
00924
00925 Double_t p2 = fP[3]*fP[3] + fP[4]*fP[4] + fP[5]*fP[5];
00926 if( p2<1.e-4 ) p2 = 1;
00927 return ( fP[3]*(xyz[0]-fP[0]) + fP[4]*(xyz[1]-fP[1]) + fP[5]*(xyz[2]-fP[2]) )/p2;
00928 }
00929
00930
00931 Double_t KFParticleBase::GetDStoPointBz( Double_t B, const Double_t xyz[] )
00932 const
00933 {
00934
00935
00936 const Double_t kCLight = 0.000299792458;
00937 Double_t bq = B*fQ*kCLight;
00938 Double_t pt2 = fP[3]*fP[3] + fP[4]*fP[4];
00939 if( pt2<1.e-4 ) return 0;
00940 Double_t dx = xyz[0] - fP[0];
00941 Double_t dy = xyz[1] - fP[1];
00942 Double_t a = dx*fP[3]+dy*fP[4];
00943 Double_t dS;
00944
00945 if( TMath::Abs(bq)<1.e-8 ) dS = a/pt2;
00946 else dS = TMath::ATan2( bq*a, pt2 + bq*(dy*fP[3] -dx*fP[4]) )/bq;
00947
00948 if(0){
00949
00950 Double_t px = fP[3];
00951 Double_t py = fP[4];
00952 Double_t pz = fP[5];
00953 Double_t ss[2], g[2][5];
00954
00955 ss[0] = dS;
00956 ss[1] = -dS;
00957 for( Int_t i=0; i<2; i++){
00958 Double_t bs = bq*ss[i];
00959 Double_t c = TMath::Cos(bs), s = TMath::Sin(bs);
00960 Double_t cB,sB;
00961 if( TMath::Abs(bq)>1.e-8){
00962 cB= (1-c)/bq;
00963 sB= s/bq;
00964 }else{
00965 const Double_t kOvSqr6 = 1./TMath::Sqrt(6.);
00966 sB = (1.-bs*kOvSqr6)*(1.+bs*kOvSqr6)*ss[i];
00967 cB = .5*sB*bs;
00968 }
00969 g[i][0] = fP[0] + sB*px + cB*py;
00970 g[i][1] = fP[1] - cB*px + sB*py;
00971 g[i][2] = fP[2] + ss[i]*pz;
00972 g[i][3] = + c*px + s*py;
00973 g[i][4] = - s*px + c*py;
00974 }
00975
00976 Int_t i=0;
00977
00978 Double_t dMin = 1.e10;
00979 for( Int_t j=0; j<2; j++){
00980 Double_t xx = g[j][0]-xyz[0];
00981 Double_t yy = g[j][1]-xyz[1];
00982 Double_t zz = g[j][2]-xyz[2];
00983 Double_t d = xx*xx + yy*yy + zz*zz;
00984 if( d<dMin ){
00985 dMin = d;
00986 i = j;
00987 }
00988 }
00989
00990 dS = ss[i];
00991
00992 Double_t x= g[i][0], y= g[i][1], z= g[i][2], ppx= g[i][3], ppy= g[i][4];
00993 Double_t ddx = x-xyz[0];
00994 Double_t ddy = y-xyz[1];
00995 Double_t ddz = z-xyz[2];
00996 Double_t c = ddx*ppx + ddy*ppy + ddz*pz ;
00997 Double_t pp2 = ppx*ppx + ppy*ppy + pz*pz;
00998 if( TMath::Abs(pp2)>1.e-8 ){
00999 dS+=c/pp2;
01000 }
01001 }
01002 return dS;
01003 }
01004
01005
01006 void KFParticleBase::GetDStoParticleBz( Double_t B, const KFParticleBase &p,
01007 Double_t &DS, Double_t &DS1 )
01008 const
01009 {
01010
01011 Double_t px = fP[3];
01012 Double_t py = fP[4];
01013 Double_t pz = fP[5];
01014
01015 Double_t px1 = p.fP[3];
01016 Double_t py1 = p.fP[4];
01017 Double_t pz1 = p.fP[5];
01018
01019 const Double_t kCLight = 0.000299792458;
01020
01021 Double_t bq = B*fQ*kCLight;
01022 Double_t bq1 = B*p.fQ*kCLight;
01023 Double_t s=0, ds=0, s1=0, ds1=0;
01024
01025 if( TMath::Abs(bq)>1.e-8 || TMath::Abs(bq1)>1.e-8 ){
01026
01027 Double_t dx = (p.fP[0] - fP[0]);
01028 Double_t dy = (p.fP[1] - fP[1]);
01029 Double_t d2 = (dx*dx+dy*dy);
01030
01031 Double_t p2 = (px *px + py *py);
01032 Double_t p21 = (px1*px1 + py1*py1);
01033
01034 Double_t a = (px*py1 - py*px1);
01035 Double_t b = (px*px1 + py*py1);
01036
01037 Double_t ldx = bq*bq1*dx - bq1*py + bq*py1 ;
01038 Double_t ldy = bq*bq1*dy + bq1*px - bq*px1 ;
01039 Double_t l2 = ldx*ldx + ldy*ldy;
01040
01041 Double_t cS = bq1*p2 + bq*bq1*(dy* px - dx* py) - bq*b;
01042 Double_t cS1= bq*p21 - bq*bq1*(dy*px1 - dx*py1) - bq1*b;
01043
01044 Double_t ca = bq*bq*bq1*d2 +2*( cS + bq*bq*(py1*dx-px1*dy)) ;
01045 Double_t ca1 = bq*bq1*bq1*d2 +2*( cS1 - bq1*bq1*(py*dx-px*dy)) ;
01046
01047 Double_t sa = 4*l2*p2 - ca*ca;
01048 Double_t sa1 = 4*l2*p21 - ca1*ca1;
01049
01050 if(sa<0) sa=0;
01051 if(sa1<0)sa1=0;
01052
01053 if( TMath::Abs(bq)>1.e-8){
01054 s = TMath::ATan2( bq*( bq1*(dx*px +dy*py) + a ) , cS )/bq;
01055 ds = TMath::ATan2(TMath::Sqrt(sa),ca)/bq;
01056 } else {
01057 s = ( (dx*px + dy*py) + (py*px1-px*py1)/bq1)/p2;
01058 ds = s*s - (d2-2*(px1*dy-py1*dx)/bq1)/p2;
01059 if( ds<0 ) ds = 0;
01060 ds = TMath::Sqrt(ds);
01061 }
01062
01063 if( TMath::Abs(bq1)>1.e-8){
01064 s1 = TMath::ATan2( -bq1*( bq*(dx*px1+dy*py1) + a), cS1 )/bq1;
01065 ds1 = TMath::ATan2(TMath::Sqrt(sa1),ca1)/bq1;
01066 } else {
01067 s1 = (-(dx*px1 + dy*py1) + (py*px1-px*py1)/bq)/p21;
01068 ds1 = s1*s1 - (d2+2*(px*dy-py*dx)/bq)/p21;
01069 if( ds1<0 ) ds1 = 0;
01070 ds1 = TMath::Sqrt(ds1);
01071 }
01072 }
01073
01074 Double_t ss[2], ss1[2], g[2][5],g1[2][5];
01075
01076 ss[0] = s + ds;
01077 ss[1] = s - ds;
01078 ss1[0] = s1 + ds1;
01079 ss1[1] = s1 - ds1;
01080 for( Int_t i=0; i<2; i++){
01081 Double_t bs = bq*ss[i];
01082 Double_t c = TMath::Cos(bs), sss = TMath::Sin(bs);
01083 Double_t cB,sB;
01084 if( TMath::Abs(bq)>1.e-8){
01085 cB= (1-c)/bq;
01086 sB= sss/bq;
01087 }else{
01088 const Double_t kOvSqr6 = 1./TMath::Sqrt(6.);
01089 sB = (1.-bs*kOvSqr6)*(1.+bs*kOvSqr6)*ss[i];
01090 cB = .5*sB*bs;
01091 }
01092 g[i][0] = fP[0] + sB*px + cB*py;
01093 g[i][1] = fP[1] - cB*px + sB*py;
01094 g[i][2] = fP[2] + ss[i]*pz;
01095 g[i][3] = + c*px + sss*py;
01096 g[i][4] = - sss*px + c*py;
01097
01098 bs = bq1*ss1[i];
01099 c = TMath::Cos(bs); sss = TMath::Sin(bs);
01100 if( TMath::Abs(bq1)>1.e-8){
01101 cB= (1-c)/bq1;
01102 sB= sss/bq1;
01103 }else{
01104 const Double_t kOvSqr6 = 1./TMath::Sqrt(6.);
01105 sB = (1.-bs*kOvSqr6)*(1.+bs*kOvSqr6)*ss1[i];
01106 cB = .5*sB*bs;
01107 }
01108
01109 g1[i][0] = p.fP[0] + sB*px1 + cB*py1;
01110 g1[i][1] = p.fP[1] - cB*px1 + sB*py1;
01111 g1[i][2] = p.fP[2] + ss[i]*pz1;
01112 g1[i][3] = + c*px1 + sss*py1;
01113 g1[i][4] = - sss*px1 + c*py1;
01114 }
01115
01116 Int_t i=0, i1=0;
01117
01118 Double_t dMin = 1.e10;
01119 for( Int_t j=0; j<2; j++){
01120 for( Int_t j1=0; j1<2; j1++){
01121 Double_t xx = g[j][0]-g1[j1][0];
01122 Double_t yy = g[j][1]-g1[j1][1];
01123 Double_t zz = g[j][2]-g1[j1][2];
01124 Double_t d = xx*xx + yy*yy + zz*zz;
01125 if( d<dMin ){
01126 dMin = d;
01127 i = j;
01128 i1 = j1;
01129 }
01130 }
01131 }
01132
01133 DS = ss[i];
01134 DS1 = ss1[i1];
01135 if(0){
01136 Double_t x= g[i][0], y= g[i][1], z= g[i][2], ppx= g[i][3], ppy= g[i][4];
01137 Double_t x1=g1[i1][0], y1= g1[i1][1], z1= g1[i1][2], ppx1= g1[i1][3], ppy1= g1[i1][4];
01138 Double_t dx = x1-x;
01139 Double_t dy = y1-y;
01140 Double_t dz = z1-z;
01141 Double_t a = ppx*ppx1 + ppy*ppy1 + pz*pz1;
01142 Double_t b = dx*ppx1 + dy*ppy1 + dz*pz1;
01143 Double_t c = dx*ppx + dy*ppy + dz*pz ;
01144 Double_t pp2 = ppx*ppx + ppy*ppy + pz*pz;
01145 Double_t pp21= ppx1*ppx1 + ppy1*ppy1 + pz1*pz1;
01146 Double_t det = pp2*pp21 - a*a;
01147 if( TMath::Abs(det)>1.e-8 ){
01148 DS+=(a*b-pp21*c)/det;
01149 DS1+=(a*c-pp2*b)/det;
01150 }
01151 }
01152 }
01153
01154
01155
01156 void KFParticleBase::TransportCBM( Double_t dS,
01157 Double_t P[], Double_t C[] ) const
01158 {
01159
01160
01161 if( fQ==0 ){
01162 TransportLine( dS, P, C );
01163 return;
01164 }
01165
01166 const Double_t kCLight = 0.000299792458;
01167
01168 Double_t c = fQ*kCLight;
01169
01170
01171
01172 Double_t
01173 px = fP[3],
01174 py = fP[4],
01175 pz = fP[5];
01176
01177 Double_t sx=0, sy=0, sz=0, syy=0, syz=0, syyy=0, ssx=0, ssy=0, ssz=0, ssyy=0, ssyz=0, ssyyy=0;
01178
01179 {
01180
01181 Double_t fld[3][3];
01182 Double_t p0[3], p1[3], p2[3];
01183
01184
01185
01186 p0[0] = fP[0];
01187 p0[1] = fP[1];
01188 p0[2] = fP[2];
01189
01190 p2[0] = fP[0] + px*dS;
01191 p2[1] = fP[1] + py*dS;
01192 p2[2] = fP[2] + pz*dS;
01193
01194 p1[0] = 0.5*(p0[0]+p2[0]);
01195 p1[1] = 0.5*(p0[1]+p2[1]);
01196 p1[2] = 0.5*(p0[2]+p2[2]);
01197
01198
01199 {
01200 GetFieldValue( p0, fld[0] );
01201 GetFieldValue( p1, fld[1] );
01202 GetFieldValue( p2, fld[2] );
01203
01204 Double_t ssy1 = ( 7*fld[0][1] + 6*fld[1][1]-fld[2][1] )*c*dS*dS/96.;
01205 Double_t ssy2 = ( fld[0][1] + 2*fld[1][1] )*c*dS*dS/6.;
01206
01207 p1[0] -= ssy1*pz;
01208 p1[2] += ssy1*px;
01209 p2[0] -= ssy2*pz;
01210 p2[2] += ssy2*px;
01211 }
01212
01213 GetFieldValue( p0, fld[0] );
01214 GetFieldValue( p1, fld[1] );
01215 GetFieldValue( p2, fld[2] );
01216
01217 sx = c*( fld[0][0] + 4*fld[1][0] + fld[2][0] )*dS/6.;
01218 sy = c*( fld[0][1] + 4*fld[1][1] + fld[2][1] )*dS/6.;
01219 sz = c*( fld[0][2] + 4*fld[1][2] + fld[2][2] )*dS/6.;
01220
01221 ssx = c*( fld[0][0] + 2*fld[1][0])*dS*dS/6.;
01222 ssy = c*( fld[0][1] + 2*fld[1][1])*dS*dS/6.;
01223 ssz = c*( fld[0][2] + 2*fld[1][2])*dS*dS/6.;
01224
01225 Double_t c2[3][3] = { { 5, -4, -1},{ 44, 80, -4},{ 11, 44, 5} };
01226 Double_t cc2[3][3] = { { 38, 8, -4},{ 148, 208, -20},{ 3, 36, 3} };
01227 for(Int_t n=0; n<3; n++)
01228 for(Int_t m=0; m<3; m++)
01229 {
01230 syz += c2[n][m]*fld[n][1]*fld[m][2];
01231 ssyz += cc2[n][m]*fld[n][1]*fld[m][2];
01232 }
01233
01234 syz *= c*c*dS*dS/360.;
01235 ssyz *= c*c*dS*dS*dS/2520.;
01236
01237 syy = c*( fld[0][1] + 4*fld[1][1] + fld[2][1] )*dS;
01238 syyy = syy*syy*syy / 1296;
01239 syy = syy*syy/72;
01240
01241 ssyy = ( fld[0][1]*( 38*fld[0][1] + 156*fld[1][1] - fld[2][1] )+
01242 fld[1][1]*( 208*fld[1][1] +16*fld[2][1] )+
01243 fld[2][1]*( 3*fld[2][1] )
01244 )*dS*dS*dS*c*c/2520.;
01245 ssyyy =
01246 (
01247 fld[0][1]*( fld[0][1]*( 85*fld[0][1] + 526*fld[1][1] - 7*fld[2][1] )+
01248 fld[1][1]*( 1376*fld[1][1] +84*fld[2][1] )+
01249 fld[2][1]*( 19*fld[2][1] ) )+
01250 fld[1][1]*( fld[1][1]*( 1376*fld[1][1] +256*fld[2][1] )+
01251 fld[2][1]*( 62*fld[2][1] ) )+
01252 fld[2][1]*fld[2][1] *( 3*fld[2][1] )
01253 )*dS*dS*dS*dS*c*c*c/90720.;
01254
01255 }
01256
01257 Double_t mJ[8][8];
01258 for( Int_t i=0; i<8; i++ ) for( Int_t j=0; j<8; j++) mJ[i][j]=0;
01259
01260 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;
01261 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;
01262 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;
01263
01264 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;
01265 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;
01266 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;
01267 mJ[6][6] = mJ[7][7] = 1;
01268
01269 P[0] = fP[0] + mJ[0][3]*px + mJ[0][4]*py + mJ[0][5]*pz;
01270 P[1] = fP[1] + mJ[1][3]*px + mJ[1][4]*py + mJ[1][5]*pz;
01271 P[2] = fP[2] + mJ[2][3]*px + mJ[2][4]*py + mJ[2][5]*pz;
01272 P[3] = mJ[3][3]*px + mJ[3][4]*py + mJ[3][5]*pz;
01273 P[4] = mJ[4][3]*px + mJ[4][4]*py + mJ[4][5]*pz;
01274 P[5] = mJ[5][3]*px + mJ[5][4]*py + mJ[5][5]*pz;
01275 P[6] = fP[6];
01276 P[7] = fP[7];
01277
01278 MultQSQt( mJ[0], fC, C);
01279
01280 }
01281
01282
01283 void KFParticleBase::TransportBz( Double_t b, Double_t t,
01284 Double_t p[], Double_t e[] ) const
01285 {
01286
01287
01288 const Double_t kCLight = 0.000299792458;
01289 b = b*fQ*kCLight;
01290 Double_t bs= b*t;
01291 Double_t s = TMath::Sin(bs), c = TMath::Cos(bs);
01292 Double_t sB, cB;
01293 if( TMath::Abs(bs)>1.e-10){
01294 sB= s/b;
01295 cB= (1-c)/b;
01296 }else{
01297 const Double_t kOvSqr6 = 1./TMath::Sqrt(6.);
01298 sB = (1.-bs*kOvSqr6)*(1.+bs*kOvSqr6)*t;
01299 cB = .5*sB*bs;
01300 }
01301
01302 Double_t px = fP[3];
01303 Double_t py = fP[4];
01304 Double_t pz = fP[5];
01305
01306 p[0] = fP[0] + sB*px + cB*py;
01307 p[1] = fP[1] - cB*px + sB*py;
01308 p[2] = fP[2] + t*pz;
01309 p[3] = c*px + s*py;
01310 p[4] = -s*px + c*py;
01311 p[5] = fP[5];
01312 p[6] = fP[6];
01313 p[7] = fP[7];
01314
01315
01316
01317
01318
01319
01320
01321
01322
01323
01324
01325
01326
01327
01328
01329
01330
01331
01332
01333
01334
01335
01336
01337
01338
01339
01340
01341
01342
01343
01344 Double_t
01345 c6=fC[6], c7=fC[7], c8=fC[8], c17=fC[17], c18=fC[18],
01346 c24 = fC[24], c31 = fC[31];
01347
01348 Double_t
01349 cBC13 = cB*fC[13],
01350 mJC13 = c7 - cB*fC[9] + sB*fC[13],
01351 mJC14 = fC[11] - cBC13 + sB*fC[14],
01352 mJC23 = c8 + t*c18,
01353 mJC24 = fC[12] + t*fC[19],
01354 mJC33 = c*fC[9] + s*fC[13],
01355 mJC34 = c*fC[13] + s*fC[14],
01356 mJC43 = -s*fC[9] + c*fC[13],
01357 mJC44 = -s*fC[13] + c*fC[14];
01358
01359
01360 e[0]= fC[0] + 2*(sB*c6 + cB*fC[10]) + (sB*fC[9] + 2*cBC13)*sB + cB*cB*fC[14];
01361 e[1]= fC[1] - cB*c6 + sB*fC[10] + mJC13*sB + mJC14*cB;
01362 e[2]= fC[2] - cB*c7 + sB*fC[11] - mJC13*cB + mJC14*sB;
01363 e[3]= fC[3] + t*fC[15] + mJC23*sB + mJC24*cB;
01364 e[4]= fC[4] + t*fC[16] - mJC23*cB + mJC24*sB;
01365
01366 e[15]= fC[15] + c18*sB + fC[19]*cB;
01367 e[16]= fC[16] - c18*cB + fC[19]*sB;
01368 e[17]= c17 + fC[20]*t;
01369 e[18]= c18*c + fC[19]*s;
01370 e[19]= -c18*s + fC[19]*c;
01371
01372 e[5]= fC[5] + (c17 + e[17] )*t;
01373
01374 e[6]= c*c6 + s*fC[10] + mJC33*sB + mJC34*cB;
01375 e[7]= c*c7 + s*fC[11] - mJC33*cB + mJC34*sB;
01376 e[8]= c*c8 + s*fC[12] + e[18]*t;
01377 e[9]= mJC33*c + mJC34*s;
01378 e[10]= -s*c6 + c*fC[10] + mJC43*sB + mJC44*cB;
01379
01380
01381 e[11]= -s*c7 + c*fC[11] - mJC43*cB + mJC44*sB;
01382 e[12]= -s*c8 + c*fC[12] + e[19]*t;
01383 e[13]= mJC43*c + mJC44*s;
01384 e[14]= -mJC43*s + mJC44*c;
01385 e[20]= fC[20];
01386 e[21]= fC[21] + fC[25]*cB + c24*sB;
01387 e[22]= fC[22] - c24*cB + fC[25]*sB;
01388 e[23]= fC[23] + fC[26]*t;
01389 e[24]= c*c24 + s*fC[25];
01390 e[25]= c*fC[25] - c24*s;
01391 e[26]= fC[26];
01392 e[27]= fC[27];
01393 e[28]= fC[28] + fC[32]*cB + c31*sB;
01394 e[29]= fC[29] - c31*cB + fC[32]*sB;
01395 e[30]= fC[30] + fC[33]*t;
01396 e[31]= c*c31 + s*fC[32];
01397 e[32]= c*fC[32] - s*c31;
01398 e[33]= fC[33];
01399 e[34]= fC[34];
01400 e[35]= fC[35];
01401 }
01402
01403
01404 Double_t KFParticleBase::GetDistanceFromVertex( const KFParticleBase &Vtx ) const
01405 {
01406
01407
01408 return GetDistanceFromVertex( Vtx.fP );
01409 }
01410
01411 Double_t KFParticleBase::GetDistanceFromVertex( const Double_t vtx[] ) const
01412 {
01413
01414
01415 Double_t mP[8], mC[36];
01416 Transport( GetDStoPoint(vtx), mP, mC );
01417 Double_t d[3]={ vtx[0]-mP[0], vtx[1]-mP[1], vtx[2]-mP[2]};
01418 return TMath::Sqrt( d[0]*d[0]+d[1]*d[1]+d[2]*d[2] );
01419 }
01420
01421 Double_t KFParticleBase::GetDistanceFromParticle( const KFParticleBase &p )
01422 const
01423 {
01424
01425
01426 Double_t dS, dS1;
01427 GetDStoParticle( p, dS, dS1 );
01428 Double_t mP[8], mC[36], mP1[8], mC1[36];
01429 Transport( dS, mP, mC );
01430 p.Transport( dS1, mP1, mC1 );
01431 Double_t dx = mP[0]-mP1[0];
01432 Double_t dy = mP[1]-mP1[1];
01433 Double_t dz = mP[2]-mP1[2];
01434 return TMath::Sqrt(dx*dx+dy*dy+dz*dz);
01435 }
01436
01437 Double_t KFParticleBase::GetDeviationFromVertex( const KFParticleBase &Vtx ) const
01438 {
01439
01440
01441 return GetDeviationFromVertex( Vtx.fP, Vtx.fC );
01442 }
01443
01444
01445 Double_t KFParticleBase::GetDeviationFromVertex( const Double_t v[], const Double_t Cv[] ) const
01446 {
01447
01448
01449
01450 Double_t mP[8];
01451 Double_t mC[36];
01452 Double_t dS = GetDStoPoint(v);
01453 Transport(dS , mP, mC );
01454
01455 Double_t d[3]={ v[0]-mP[0], v[1]-mP[1], v[2]-mP[2]};
01456
01457 Double_t sigmaS = .1+10.*TMath::Sqrt( (d[0]*d[0]+d[1]*d[1]+d[2]*d[2])/
01458 (mP[3]*mP[3]+mP[4]*mP[4]+mP[5]*mP[5]) );
01459
01460
01461 Double_t h[3] = { mP[3]*sigmaS, mP[4]*sigmaS, mP[5]*sigmaS };
01462
01463 Double_t mSi[6] =
01464 { mC[0] +h[0]*h[0],
01465 mC[1] +h[1]*h[0], mC[2] +h[1]*h[1],
01466 mC[3] +h[2]*h[0], mC[4] +h[2]*h[1], mC[5] +h[2]*h[2] };
01467
01468 if( Cv ){
01469 mSi[0]+=Cv[0];
01470 mSi[1]+=Cv[1];
01471 mSi[2]+=Cv[2];
01472 mSi[3]+=Cv[3];
01473 mSi[4]+=Cv[4];
01474 mSi[5]+=Cv[5];
01475 }
01476
01477 Double_t mS[6];
01478
01479 mS[0] = mSi[2]*mSi[5] - mSi[4]*mSi[4];
01480 mS[1] = mSi[3]*mSi[4] - mSi[1]*mSi[5];
01481 mS[2] = mSi[0]*mSi[5] - mSi[3]*mSi[3];
01482 mS[3] = mSi[1]*mSi[4] - mSi[2]*mSi[3];
01483 mS[4] = mSi[1]*mSi[3] - mSi[0]*mSi[4];
01484 mS[5] = mSi[0]*mSi[2] - mSi[1]*mSi[1];
01485
01486 Double_t s = ( mSi[0]*mS[0] + mSi[1]*mS[1] + mSi[3]*mS[3] );
01487 s = ( s > 1.E-20 ) ?1./s :0;
01488
01489 return TMath::Sqrt( TMath::Abs(s*( ( mS[0]*d[0] + mS[1]*d[1] + mS[3]*d[2])*d[0]
01490 +(mS[1]*d[0] + mS[2]*d[1] + mS[4]*d[2])*d[1]
01491 +(mS[3]*d[0] + mS[4]*d[1] + mS[5]*d[2])*d[2] ))/2);
01492 }
01493
01494
01495 Double_t KFParticleBase::GetDeviationFromParticle( const KFParticleBase &p )
01496 const
01497 {
01498
01499
01500 Double_t dS, dS1;
01501 GetDStoParticle( p, dS, dS1 );
01502 Double_t mP1[8], mC1[36];
01503 p.Transport( dS1, mP1, mC1 );
01504
01505 Double_t d[3]={ fP[0]-mP1[0], fP[1]-mP1[1], fP[2]-mP1[2]};
01506
01507 Double_t sigmaS = .1+10.*TMath::Sqrt( (d[0]*d[0]+d[1]*d[1]+d[2]*d[2])/
01508 (mP1[3]*mP1[3]+mP1[4]*mP1[4]+mP1[5]*mP1[5]) );
01509
01510 Double_t h[3] = { mP1[3]*sigmaS, mP1[4]*sigmaS, mP1[5]*sigmaS };
01511
01512 mC1[0] +=h[0]*h[0];
01513 mC1[1] +=h[1]*h[0];
01514 mC1[2] +=h[1]*h[1];
01515 mC1[3] +=h[2]*h[0];
01516 mC1[4] +=h[2]*h[1];
01517 mC1[5] +=h[2]*h[2];
01518
01519 return GetDeviationFromVertex( mP1, mC1 )*TMath::Sqrt(2./1.);
01520 }
01521
01522
01523
01524 void KFParticleBase::SubtractFromVertex( KFParticleBase &Vtx ) const
01525 {
01526
01527
01528 Double_t fld[3];
01529 {
01530 GetFieldValue( Vtx.fP, fld );
01531 const Double_t kCLight = 0.000299792458;
01532 fld[0]*=kCLight; fld[1]*=kCLight; fld[2]*=kCLight;
01533 }
01534
01535 Double_t m[8];
01536 Double_t mCm[36];
01537
01538 if( Vtx.fIsLinearized ){
01539 GetMeasurement( Vtx.fVtxGuess, m, mCm );
01540 } else {
01541 GetMeasurement( Vtx.fP, m, mCm );
01542 }
01543
01544 Double_t mV[6];
01545
01546 mV[ 0] = mCm[ 0];
01547 mV[ 1] = mCm[ 1];
01548 mV[ 2] = mCm[ 2];
01549 mV[ 3] = mCm[ 3];
01550 mV[ 4] = mCm[ 4];
01551 mV[ 5] = mCm[ 5];
01552
01553
01554
01555 Double_t mS[6];
01556 {
01557 Double_t mSi[6] = { mV[0]-Vtx.fC[0],
01558 mV[1]-Vtx.fC[1], mV[2]-Vtx.fC[2],
01559 mV[3]-Vtx.fC[3], mV[4]-Vtx.fC[4], mV[5]-Vtx.fC[5] };
01560
01561 mS[0] = mSi[2]*mSi[5] - mSi[4]*mSi[4];
01562 mS[1] = mSi[3]*mSi[4] - mSi[1]*mSi[5];
01563 mS[2] = mSi[0]*mSi[5] - mSi[3]*mSi[3];
01564 mS[3] = mSi[1]*mSi[4] - mSi[2]*mSi[3];
01565 mS[4] = mSi[1]*mSi[3] - mSi[0]*mSi[4];
01566 mS[5] = mSi[0]*mSi[2] - mSi[1]*mSi[1];
01567
01568 Double_t s = ( mSi[0]*mS[0] + mSi[1]*mS[1] + mSi[3]*mS[3] );
01569 s = ( s > 1.E-20 ) ?1./s :0;
01570 mS[0]*=s;
01571 mS[1]*=s;
01572 mS[2]*=s;
01573 mS[3]*=s;
01574 mS[4]*=s;
01575 mS[5]*=s;
01576 }
01577
01578
01579
01580 Double_t zeta[3] = { m[0]-Vtx.fP[0], m[1]-Vtx.fP[1], m[2]-Vtx.fP[2] };
01581
01582
01583
01584 Double_t mCHt0[3], mCHt1[3], mCHt2[3];
01585
01586 mCHt0[0]=Vtx.fC[ 0] ; mCHt1[0]=Vtx.fC[ 1] ; mCHt2[0]=Vtx.fC[ 3] ;
01587 mCHt0[1]=Vtx.fC[ 1] ; mCHt1[1]=Vtx.fC[ 2] ; mCHt2[1]=Vtx.fC[ 4] ;
01588 mCHt0[2]=Vtx.fC[ 3] ; mCHt1[2]=Vtx.fC[ 4] ; mCHt2[2]=Vtx.fC[ 5] ;
01589
01590
01591
01592 Double_t k0[3], k1[3], k2[3];
01593
01594 for(Int_t i=0;i<3;++i){
01595 k0[i] = mCHt0[i]*mS[0] + mCHt1[i]*mS[1] + mCHt2[i]*mS[3];
01596 k1[i] = mCHt0[i]*mS[1] + mCHt1[i]*mS[2] + mCHt2[i]*mS[4];
01597 k2[i] = mCHt0[i]*mS[3] + mCHt1[i]*mS[4] + mCHt2[i]*mS[5];
01598 }
01599
01600
01601
01602 Double_t dChi2 = -(mS[0]*zeta[0] + mS[1]*zeta[1] + mS[3]*zeta[2])*zeta[0]
01603 + (mS[1]*zeta[0] + mS[2]*zeta[1] + mS[4]*zeta[2])*zeta[1]
01604 + (mS[3]*zeta[0] + mS[4]*zeta[1] + mS[5]*zeta[2])*zeta[2];
01605
01606 if( Vtx.fChi2 - dChi2 < 0 ) return;
01607
01608 for(Int_t i=0;i<3;++i)
01609 Vtx.fP[i] -= k0[i]*zeta[0] + k1[i]*zeta[1] + k2[i]*zeta[2];
01610
01611
01612
01613 for(Int_t i=0, k=0;i<3;++i){
01614 for(Int_t j=0;j<=i;++j,++k)
01615 Vtx.fC[k] += k0[i]*mCHt0[j] + k1[i]*mCHt1[j] + k2[i]*mCHt2[j];
01616 }
01617
01618
01619
01620 Vtx.fNDF -= 2;
01621 Vtx.fChi2 -= dChi2;
01622 }
01623
01624
01625
01626 void KFParticleBase::TransportLine( Double_t dS,
01627 Double_t P[], Double_t C[] ) const
01628 {
01629
01630
01631 P[0] = fP[0] + dS*fP[3];
01632 P[1] = fP[1] + dS*fP[4];
01633 P[2] = fP[2] + dS*fP[5];
01634 P[3] = fP[3];
01635 P[4] = fP[4];
01636 P[5] = fP[5];
01637 P[6] = fP[6];
01638 P[7] = fP[7];
01639
01640 Double_t c6 = fC[ 6] + dS*fC[ 9];
01641 Double_t c11 = fC[11] + dS*fC[14];
01642 Double_t c17 = fC[17] + dS*fC[20];
01643 Double_t sc13 = dS*fC[13];
01644 Double_t sc18 = dS*fC[18];
01645 Double_t sc19 = dS*fC[19];
01646
01647 C[ 0] = fC[ 0] + dS*( fC[ 6] + c6 );
01648 C[ 2] = fC[ 2] + dS*( fC[11] + c11 );
01649 C[ 5] = fC[ 5] + dS*( fC[17] + c17 );
01650
01651 C[ 7] = fC[ 7] + sc13;
01652 C[ 8] = fC[ 8] + sc18;
01653 C[ 9] = fC[ 9];
01654
01655 C[12] = fC[12] + sc19;
01656
01657 C[ 1] = fC[ 1] + dS*( fC[10] + C[ 7] );
01658 C[ 3] = fC[ 3] + dS*( fC[15] + C[ 8] );
01659 C[ 4] = fC[ 4] + dS*( fC[16] + C[12] );
01660 C[ 6] = c6;
01661
01662 C[10] = fC[10] + sc13;
01663 C[11] = c11;
01664
01665 C[13] = fC[13];
01666 C[14] = fC[14];
01667 C[15] = fC[15] + sc18;
01668 C[16] = fC[16] + sc19;
01669 C[17] = c17;
01670
01671 C[18] = fC[18];
01672 C[19] = fC[19];
01673 C[20] = fC[20];
01674 C[21] = fC[21] + dS*fC[24];
01675 C[22] = fC[22] + dS*fC[25];
01676 C[23] = fC[23] + dS*fC[26];
01677
01678 C[24] = fC[24];
01679 C[25] = fC[25];
01680 C[26] = fC[26];
01681 C[27] = fC[27];
01682 C[28] = fC[28] + dS*fC[31];
01683 C[29] = fC[29] + dS*fC[32];
01684 C[30] = fC[30] + dS*fC[33];
01685
01686 C[31] = fC[31];
01687 C[32] = fC[32];
01688 C[33] = fC[33];
01689 C[34] = fC[34];
01690 C[35] = fC[35];
01691 }
01692
01693
01694 void KFParticleBase::ConstructGammaBz( const KFParticleBase &daughter1,
01695 const KFParticleBase &daughter2, Double_t Bz )
01696 {
01697
01698
01699 const KFParticleBase *daughters[2] = { &daughter1, &daughter2};
01700
01701 Double_t v0[3];
01702
01703 if( !fIsLinearized ){
01704 Double_t ds, ds1;
01705 Double_t m[8];
01706 Double_t mCd[36];
01707 daughter1.GetDStoParticle(daughter2, ds, ds1);
01708 daughter1.Transport( ds, m, mCd );
01709 v0[0] = m[0];
01710 v0[1] = m[1];
01711 v0[2] = m[2];
01712 daughter2.Transport( ds1, m, mCd );
01713 v0[0] = .5*( v0[0] + m[0] );
01714 v0[1] = .5*( v0[1] + m[1] );
01715 v0[2] = .5*( v0[2] + m[2] );
01716 } else {
01717 v0[0] = fVtxGuess[0];
01718 v0[1] = fVtxGuess[1];
01719 v0[2] = fVtxGuess[2];
01720 }
01721
01722 fAtProductionVertex = 0;
01723 fSFromDecay = 0;
01724 fP[0] = v0[0];
01725 fP[1] = v0[1];
01726 fP[2] = v0[2];
01727 fP[3] = 0;
01728 fP[4] = 0;
01729 fP[5] = 0;
01730 fP[6] = 0;
01731 fP[7] = 0;
01732
01733
01734
01735
01736 Double_t daughterP[2][8], daughterC[2][36];
01737 Double_t vtxMom[2][3];
01738
01739 {
01740 for( Int_t id=0; id<2; id++ ){
01741
01742 Double_t *p = daughterP[id];
01743 Double_t *mC = daughterC[id];
01744
01745 daughters[id]->GetMeasurement( v0, p, mC );
01746
01747 Double_t mAi[6];
01748 InvertSym3(mC, mAi );
01749
01750 Double_t mB[3][3];
01751
01752 mB[0][0] = mC[ 6]*mAi[0] + mC[ 7]*mAi[1] + mC[ 8]*mAi[3];
01753 mB[0][1] = mC[ 6]*mAi[1] + mC[ 7]*mAi[2] + mC[ 8]*mAi[4];
01754 mB[0][2] = mC[ 6]*mAi[3] + mC[ 7]*mAi[4] + mC[ 8]*mAi[5];
01755
01756 mB[1][0] = mC[10]*mAi[0] + mC[11]*mAi[1] + mC[12]*mAi[3];
01757 mB[1][1] = mC[10]*mAi[1] + mC[11]*mAi[2] + mC[12]*mAi[4];
01758 mB[1][2] = mC[10]*mAi[3] + mC[11]*mAi[4] + mC[12]*mAi[5];
01759
01760 mB[2][0] = mC[15]*mAi[0] + mC[16]*mAi[1] + mC[17]*mAi[3];
01761 mB[2][1] = mC[15]*mAi[1] + mC[16]*mAi[2] + mC[17]*mAi[4];
01762 mB[2][2] = mC[15]*mAi[3] + mC[16]*mAi[4] + mC[17]*mAi[5];
01763
01764 Double_t z[3] = { v0[0]-p[0], v0[1]-p[1], v0[2]-p[2] };
01765
01766 vtxMom[id][0] = p[3] + mB[0][0]*z[0] + mB[0][1]*z[1] + mB[0][2]*z[2];
01767 vtxMom[id][1] = p[4] + mB[1][0]*z[0] + mB[1][1]*z[1] + mB[1][2]*z[2];
01768 vtxMom[id][2] = p[5] + mB[2][0]*z[0] + mB[2][1]*z[1] + mB[2][2]*z[2];
01769
01770 daughters[id]->Transport( daughters[id]->GetDStoPoint(v0), p, mC );
01771
01772 }
01773
01774 }
01775
01776
01777
01778 {
01779
01780 Double_t mpx0 = vtxMom[0][0]+vtxMom[1][0];
01781 Double_t mpy0 = vtxMom[0][1]+vtxMom[1][1];
01782 Double_t mpt0 = TMath::Sqrt(mpx0*mpx0 + mpy0*mpy0);
01783
01784
01785 Double_t ca0 = mpx0/mpt0;
01786 Double_t sa0 = mpy0/mpt0;
01787 Double_t r[3] = { v0[0], v0[1], v0[2] };
01788 Double_t mC[3][3] = {{10000., 0 , 0 },
01789 {0, 10000., 0 },
01790 {0, 0, 10000. } };
01791 Double_t chi2=0;
01792
01793 for( Int_t id=0; id<2; id++ ){
01794 const Double_t kCLight = 0.000299792458;
01795 Double_t q = Bz*daughters[id]->GetQ()*kCLight;
01796 Double_t px0 = vtxMom[id][0];
01797 Double_t py0 = vtxMom[id][1];
01798 Double_t pz0 = vtxMom[id][2];
01799 Double_t pt0 = TMath::Sqrt(px0*px0+py0*py0);
01800 Double_t mG[3][6], mB[3], mH[3][3];
01801
01802
01803
01804
01805
01806
01807
01808
01809 mG[0][0] = q;
01810 mG[0][1] = 0;
01811 mG[0][2] = 0;
01812 mG[0][3] = -sa0*px0/pt0;
01813 mG[0][4] = 1 -sa0*py0/pt0;
01814 mG[0][5] = 0;
01815 mH[0][0] = q;
01816 mH[0][1] = 0;
01817 mH[0][2] = 0;
01818 mB[0] = py0 - sa0*pt0 - mG[0][3]*px0 - mG[0][4]*py0 ;
01819
01820
01821
01822 mG[1][0] = 0;
01823 mG[1][1] = q;
01824 mG[1][2] = 0;
01825 mG[1][3] = -1 + ca0*px0/pt0;
01826 mG[1][4] = + ca0*py0/pt0;
01827 mG[1][5] = 0;
01828 mH[1][0] = 0;
01829 mH[1][1] = q;
01830 mH[1][2] = 0;
01831 mB[1] = -px0 + ca0*pt0 - mG[1][3]*px0 - mG[1][4]*py0 ;
01832
01833
01834
01835 mG[2][0] = -pz0*ca0;
01836 mG[2][1] = -pz0*sa0;
01837 mG[2][2] = px0*ca0 + py0*sa0;
01838 mG[2][3] = 0;
01839 mG[2][4] = 0;
01840 mG[2][5] = 0;
01841
01842 mH[2][0] = mG[2][0];
01843 mH[2][1] = mG[2][1];
01844 mH[2][2] = mG[2][2];
01845
01846 mB[2] = 0;
01847
01848
01849
01850
01851
01852 Double_t mGV[3][6];
01853 Double_t mV[6];
01854 Double_t m[3];
01855 for( Int_t i=0; i<3; i++ ){
01856 m[i] = mB[i];
01857 for( Int_t k=0; k<6; k++ ) m[i]+=mG[i][k]*daughterP[id][k];
01858 }
01859 for( Int_t i=0; i<3; i++ ){
01860 for( Int_t j=0; j<6; j++ ){
01861 mGV[i][j] = 0;
01862 for( Int_t k=0; k<6; k++ ) mGV[i][j]+=mG[i][k]*daughterC[id][ IJ(k,j) ];
01863 }
01864 }
01865 for( Int_t i=0, k=0; i<3; i++ ){
01866 for( Int_t j=0; j<=i; j++,k++ ){
01867 mV[k] = 0;
01868 for( Int_t l=0; l<6; l++ ) mV[k]+=mGV[i][l]*mG[j][l];
01869 }
01870 }
01871
01872
01873
01874
01875 Double_t mCHt[3][3];
01876 Double_t mHCHt[6];
01877 Double_t mHr[3];
01878 for( Int_t i=0; i<3; i++ ){
01879 mHr[i] = 0;
01880 for( Int_t k=0; k<3; k++ ) mHr[i]+= mH[i][k]*r[k];
01881 }
01882
01883 for( Int_t i=0; i<3; i++ ){
01884 for( Int_t j=0; j<3; j++){
01885 mCHt[i][j] = 0;
01886 for( Int_t k=0; k<3; k++ ) mCHt[i][j]+= mC[i][k]*mH[j][k];
01887 }
01888 }
01889
01890 for( Int_t i=0, k=0; i<3; i++ ){
01891 for( Int_t j=0; j<=i; j++, k++ ){
01892 mHCHt[k] = 0;
01893 for( Int_t l=0; l<3; l++ ) mHCHt[k]+= mH[i][l]*mCHt[l][j];
01894 }
01895 }
01896
01897 Double_t mS[6] = { mHCHt[0]+mV[0],
01898 mHCHt[1]+mV[1], mHCHt[2]+mV[2],
01899 mHCHt[3]+mV[3], mHCHt[4]+mV[4], mHCHt[5]+mV[5] };
01900
01901
01902 InvertSym3(mS,mS);
01903
01904
01905
01906 Double_t zeta[3] = { m[0]-mHr[0], m[1]-mHr[1], m[2]-mHr[2] };
01907
01908
01909
01910 Double_t k[3][3];
01911
01912 for(Int_t i=0;i<3;++i){
01913 k[i][0] = mCHt[i][0]*mS[0] + mCHt[i][1]*mS[1] + mCHt[i][2]*mS[3];
01914 k[i][1] = mCHt[i][0]*mS[1] + mCHt[i][1]*mS[2] + mCHt[i][2]*mS[4];
01915 k[i][2] = mCHt[i][0]*mS[3] + mCHt[i][1]*mS[4] + mCHt[i][2]*mS[5];
01916 }
01917
01918
01919
01920 for(Int_t i=0;i<3;++i)
01921 r[i] = r[i] + k[i][0]*zeta[0] + k[i][1]*zeta[1] + k[i][2]*zeta[2];
01922
01923
01924
01925 for(Int_t i=0;i<3;++i){
01926 for(Int_t j=0;j<=i;++j){
01927 mC[i][j] = mC[i][j] - (k[i][0]*mCHt[j][0] + k[i][1]*mCHt[j][1] + k[i][2]*mCHt[j][2]);
01928 mC[j][i] = mC[i][j];
01929 }
01930 }
01931
01932
01933
01934
01935 chi2 += ( ( mS[0]*zeta[0] + mS[1]*zeta[1] + mS[3]*zeta[2] )*zeta[0]
01936 +(mS[1]*zeta[0] + mS[2]*zeta[1] + mS[4]*zeta[2] )*zeta[1]
01937 +(mS[3]*zeta[0] + mS[4]*zeta[1] + mS[5]*zeta[2] )*zeta[2] );
01938 }
01939
01940
01941
01942 fNDF = 2;
01943 fChi2 = chi2;
01944 for( Int_t i=0; i<3; i++ ) fP[i] = r[i];
01945 for( Int_t i=0,k=0; i<3; i++ ){
01946 for( Int_t j=0; j<=i; j++,k++ ){
01947 fC[k] = mC[i][j];
01948 }
01949 }
01950 }
01951
01952
01953
01954 fQ = 0;
01955 fSFromDecay = 0;
01956
01957 for(Int_t i=3;i<8;++i) fP[i]=0.;
01958 for(Int_t i=6;i<36;++i) fC[i]=0.;
01959
01960
01961 for( Int_t id=0; id<2; id++ ){
01962
01963 Double_t *p = daughterP[id];
01964 Double_t *mC = daughterC[id];
01965 daughters[id]->GetMeasurement( v0, p, mC );
01966
01967 const Double_t *m = fP, *mV = fC;
01968
01969 Double_t mAi[6];
01970 InvertSym3(mC, mAi );
01971
01972 Double_t mB[4][3];
01973
01974 mB[0][0] = mC[ 6]*mAi[0] + mC[ 7]*mAi[1] + mC[ 8]*mAi[3];
01975 mB[0][1] = mC[ 6]*mAi[1] + mC[ 7]*mAi[2] + mC[ 8]*mAi[4];
01976 mB[0][2] = mC[ 6]*mAi[3] + mC[ 7]*mAi[4] + mC[ 8]*mAi[5];
01977
01978 mB[1][0] = mC[10]*mAi[0] + mC[11]*mAi[1] + mC[12]*mAi[3];
01979 mB[1][1] = mC[10]*mAi[1] + mC[11]*mAi[2] + mC[12]*mAi[4];
01980 mB[1][2] = mC[10]*mAi[3] + mC[11]*mAi[4] + mC[12]*mAi[5];
01981
01982 mB[2][0] = mC[15]*mAi[0] + mC[16]*mAi[1] + mC[17]*mAi[3];
01983 mB[2][1] = mC[15]*mAi[1] + mC[16]*mAi[2] + mC[17]*mAi[4];
01984 mB[2][2] = mC[15]*mAi[3] + mC[16]*mAi[4] + mC[17]*mAi[5];
01985
01986 mB[3][0] = mC[21]*mAi[0] + mC[22]*mAi[1] + mC[23]*mAi[3];
01987 mB[3][1] = mC[21]*mAi[1] + mC[22]*mAi[2] + mC[23]*mAi[4];
01988 mB[3][2] = mC[21]*mAi[3] + mC[22]*mAi[4] + mC[23]*mAi[5];
01989
01990
01991 Double_t z[3] = { m[0]-p[0], m[1]-p[1], m[2]-p[2] };
01992
01993 {
01994 Double_t mAV[6] = { mC[0]-mV[0], mC[1]-mV[1], mC[2]-mV[2],
01995 mC[3]-mV[3], mC[4]-mV[4], mC[5]-mV[5] };
01996
01997 Double_t mAVi[6];
01998 if( !InvertSym3(mAV, mAVi) ){
01999 Double_t dChi2 = ( +(mAVi[0]*z[0] + mAVi[1]*z[1] + mAVi[3]*z[2])*z[0]
02000 +(mAVi[1]*z[0] + mAVi[2]*z[1] + mAVi[4]*z[2])*z[1]
02001 +(mAVi[3]*z[0] + mAVi[4]*z[1] + mAVi[5]*z[2])*z[2] );
02002 fChi2+= TMath::Abs( dChi2 );
02003 }
02004 fNDF += 2;
02005 }
02006
02007
02008
02009 fP[3]+= p[3] + mB[0][0]*z[0] + mB[0][1]*z[1] + mB[0][2]*z[2];
02010 fP[4]+= p[4] + mB[1][0]*z[0] + mB[1][1]*z[1] + mB[1][2]*z[2];
02011 fP[5]+= p[5] + mB[2][0]*z[0] + mB[2][1]*z[1] + mB[2][2]*z[2];
02012 fP[6]+= p[6] + mB[3][0]*z[0] + mB[3][1]*z[1] + mB[3][2]*z[2];
02013
02014 Double_t d0, d1, d2;
02015
02016 d0= mB[0][0]*mV[0] + mB[0][1]*mV[1] + mB[0][2]*mV[3] - mC[ 6];
02017 d1= mB[0][0]*mV[1] + mB[0][1]*mV[2] + mB[0][2]*mV[4] - mC[ 7];
02018 d2= mB[0][0]*mV[3] + mB[0][1]*mV[4] + mB[0][2]*mV[5] - mC[ 8];
02019
02020
02021
02022
02023 fC[9]+= mC[ 9] + d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
02024
02025 d0= mB[1][0]*mV[0] + mB[1][1]*mV[1] + mB[1][2]*mV[3] - mC[10];
02026 d1= mB[1][0]*mV[1] + mB[1][1]*mV[2] + mB[1][2]*mV[4] - mC[11];
02027 d2= mB[1][0]*mV[3] + mB[1][1]*mV[4] + mB[1][2]*mV[5] - mC[12];
02028
02029
02030
02031
02032 fC[13]+= mC[13]+ d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
02033 fC[14]+= mC[14]+ d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
02034
02035 d0= mB[2][0]*mV[0] + mB[2][1]*mV[1] + mB[2][2]*mV[3] - mC[15];
02036 d1= mB[2][0]*mV[1] + mB[2][1]*mV[2] + mB[2][2]*mV[4] - mC[16];
02037 d2= mB[2][0]*mV[3] + mB[2][1]*mV[4] + mB[2][2]*mV[5] - mC[17];
02038
02039
02040
02041
02042 fC[18]+= mC[18]+ d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
02043 fC[19]+= mC[19]+ d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
02044 fC[20]+= mC[20]+ d0*mB[2][0] + d1*mB[2][1] + d2*mB[2][2];
02045
02046 d0= mB[3][0]*mV[0] + mB[3][1]*mV[1] + mB[3][2]*mV[3] - mC[21];
02047 d1= mB[3][0]*mV[1] + mB[3][1]*mV[2] + mB[3][2]*mV[4] - mC[22];
02048 d2= mB[3][0]*mV[3] + mB[3][1]*mV[4] + mB[3][2]*mV[5] - mC[23];
02049
02050
02051
02052
02053 fC[24]+= mC[24] + d0*mB[0][0] + d1*mB[0][1] + d2*mB[0][2];
02054 fC[25]+= mC[25] + d0*mB[1][0] + d1*mB[1][1] + d2*mB[1][2];
02055 fC[26]+= mC[26] + d0*mB[2][0] + d1*mB[2][1] + d2*mB[2][2];
02056 fC[27]+= mC[27] + d0*mB[3][0] + d1*mB[3][1] + d2*mB[3][2];
02057 }
02058
02059 SetMassConstraint(0,0);
02060
02061 }
02062
02063 Bool_t KFParticleBase::InvertSym3( const Double_t A[], Double_t Ai[] )
02064 {
02065
02066
02067 Bool_t ret = 0;
02068 Double_t a0 = A[0], a1 = A[1], a2 = A[2], a3 = A[3];
02069
02070 Ai[0] = a2*A[5] - A[4]*A[4];
02071 Ai[1] = a3*A[4] - a1*A[5];
02072 Ai[3] = a1*A[4] - a2*a3;
02073 Double_t det = (a0*Ai[0] + a1*Ai[1] + a3*Ai[3]);
02074 if( det>1.e-20 ) det = 1./det;
02075 else{
02076 det = 0;
02077 ret = 1;
02078 }
02079 Ai[0] *= det;
02080 Ai[1] *= det;
02081 Ai[3] *= det;
02082 Ai[2] = ( a0*A[5] - a3*a3 )*det;
02083 Ai[4] = ( a1*a3 - a0*A[4] )*det;
02084 Ai[5] = ( a0*a2 - a1*a1 )*det;
02085 return ret;
02086 }
02087
02088 void KFParticleBase::MultQSQt( const Double_t Q[], const Double_t S[], Double_t SOut[] )
02089 {
02090
02091
02092 const Int_t kN= 8;
02093 Double_t mA[kN*kN];
02094
02095 for( Int_t i=0, ij=0; i<kN; i++ ){
02096 for( Int_t j=0; j<kN; j++, ++ij ){
02097 mA[ij] = 0 ;
02098 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];
02099 }
02100 }
02101
02102 for( Int_t i=0; i<kN; i++ ){
02103 for( Int_t j=0; j<=i; j++ ){
02104 Int_t ij = ( j<=i ) ? i*(i+1)/2+j :j*(j+1)/2+i;
02105 SOut[ij] = 0 ;
02106 for( Int_t k=0; k<kN; k++ ) SOut[ij] += Q[ i*kN+k ] * mA[ k*kN+j ];
02107 }
02108 }
02109 }
02110
02111 void KFParticleBase::Print(Option_t *opt) const {
02112 cout << *this << endl;
02113 }
02114
02115 ostream& operator<<(ostream& os, const KFParticleBase& particle) {
02116 static const Char_t *vn[11] = {"x","y","z","px","py","pz","E","S","Q","Chi2","NDF"};
02117 os << "KFP:";
02118 for (Int_t i = 0; i < 8; i++)
02119 os << Form(" %s: %7.2f +/- %7.2f", vn[i],
02120 particle.GetParameter(i),
02121 (particle.GetCovariance(i,i) >= 0) ? TMath::Sqrt(particle.GetCovariance(i,i)) : -13);
02122 os << " Q:" << particle.GetQ() << " chi2/NDF : " << particle.GetChi2() << "/" << particle.GetNDF();
02123 return os;
02124 }
02125
02126
02127
02128
02129