X-Git-Url: http://git.uio.no/git/?a=blobdiff_plain;f=ITS%2FAliITStrackV2.cxx;h=a4d86770e9bf142c443f69aa886ca6869fda3ff4;hb=20f79d37d89dc1fa42a87fab4697de017cd08cf0;hp=cf2a7f0521e16770168416a0025838e169f8949a;hpb=d3dbf32381b6ccb54730c567a989c90cd5c76bc9;p=u%2Fmrichter%2FAliRoot.git diff --git a/ITS/AliITStrackV2.cxx b/ITS/AliITStrackV2.cxx index cf2a7f0521e..a4d86770e9b 100644 --- a/ITS/AliITStrackV2.cxx +++ b/ITS/AliITStrackV2.cxx @@ -28,6 +28,7 @@ #include "AliITSReconstructor.h" #include "AliITStrackV2.h" #include "AliTracker.h" +#include "AliLog.h" const Int_t AliITStrackV2::fgkWARN = 5; @@ -41,12 +42,13 @@ AliITStrackV2::AliITStrackV2() : AliKalmanTrack(), fESDtrack(0) { for(Int_t i=0; i<2*AliITSgeomTGeo::kNLayers; i++) {fIndex[i]=-1; fModule[i]=-1;} + for(Int_t i=0; iGetX(),par->GetAlpha(),par->GetParameter(),par->GetCovariance()); - //if (!Invariant()) throw "AliITStrackV2: conversion failed !\n"; - SetLabel(t.GetLabel()); SetMass(t.GetMass()); SetNumberOfClusters(t.GetITSclusters(fIndex)); @@ -75,6 +75,7 @@ AliITStrackV2::AliITStrackV2(AliESDtrack& t,Bool_t c) throw (const Char_t *) : SetIntegratedLength(t.GetIntegratedLength()); } + for(Int_t i=0; iUpdateTrackParams(this,flags); // copy the module indices - for(Int_t i=0;i<12;i++) { + Int_t i; + for(i=0;i<2*AliITSgeomTGeo::kNLayers;i++) { // printf(" %d\n",GetModuleIndex(i)); fESDtrack->SetITSModuleIndex(i,GetModuleIndex(i)); } + // copy the map of shared clusters + if(flags==AliESDtrack::kITSin) { + UChar_t itsSharedMap=0; + for(i=0;i0) SETBIT(itsSharedMap,i); + + } + fESDtrack->SetITSSharedMap(itsSharedMap); + } + // copy the 4 dedx samples Double_t sdedx[4]={0.,0.,0.,0.}; - for(Int_t i=0; i<4; i++) sdedx[i]=fdEdxSample[i]; + for(i=0; i<4; i++) sdedx[i]=fdEdxSample[i]; fESDtrack->SetITSdEdxSamples(sdedx); } @@ -119,6 +132,7 @@ AliITStrackV2::AliITStrackV2(const AliITStrackV2& t) : fIndex[i]=t.fIndex[i]; fModule[i]=t.fModule[i]; } + for (i=0; ifgkWARN) AliWarning("Wrong invariant !"); + if (n>fgkWARN) AliDebug(1,"Wrong invariant !"); return kFALSE; } @@ -286,12 +300,16 @@ Bool_t AliITStrackV2::Invariant() const { if(!fCheckInvariant) return kTRUE; Int_t n=GetNumberOfClusters(); - + static Float_t bz = GetBz(); // take into account the misalignment error Float_t maxMisalErrY2=0,maxMisalErrZ2=0; + //RS + const AliITSRecoParam* recopar = AliITSReconstructor::GetRecoParam(); + if (!recopar) recopar = AliITSRecoParam::GetHighFluxParam(); + for (Int_t lay=0; layGetClusterMisalErrorY(lay,GetBz())); - maxMisalErrZ2 = TMath::Max(maxMisalErrZ2,AliITSReconstructor::GetRecoParam()->GetClusterMisalErrorZ(lay,GetBz())); + maxMisalErrY2 = TMath::Max(maxMisalErrY2,recopar->GetClusterMisalErrorY(lay,bz)); + maxMisalErrZ2 = TMath::Max(maxMisalErrZ2,recopar->GetClusterMisalErrorZ(lay,bz)); } maxMisalErrY2 *= maxMisalErrY2; maxMisalErrZ2 *= maxMisalErrZ2; @@ -302,32 +320,32 @@ Bool_t AliITStrackV2::Invariant() const { Double_t sP2=GetParameter()[2]; if (TMath::Abs(sP2) >= kAlmost1){ - if (n>fgkWARN) Warning("Invariant","fP2=%f\n",sP2); + if (n>fgkWARN) AliDebug(1,Form("fP2=%f\n",sP2)); return kFALSE; } Double_t sC00=GetCovariance()[0]; if (sC00<=0 || sC00>(9.+maxMisalErrY2)) { - if (n>fgkWARN) Warning("Invariant","fC00=%f\n",sC00); + if (n>fgkWARN) AliDebug(1,Form("fC00=%f\n",sC00)); return kFALSE; } Double_t sC11=GetCovariance()[2]; if (sC11<=0 || sC11>(9.+maxMisalErrZ2)) { - if (n>fgkWARN) Warning("Invariant","fC11=%f\n",sC11); + if (n>fgkWARN) AliDebug(1,Form("fC11=%f\n",sC11)); return kFALSE; } Double_t sC22=GetCovariance()[5]; if (sC22<=0 || sC22>1.) { - if (n>fgkWARN) Warning("Invariant","fC22=%f\n",sC22); + if (n>fgkWARN) AliDebug(1,Form("fC22=%f\n",sC22)); return kFALSE; } Double_t sC33=GetCovariance()[9]; if (sC33<=0 || sC33>1.) { - if (n>fgkWARN) Warning("Invariant","fC33=%f\n",sC33); + if (n>fgkWARN) AliDebug(1,Form("fC33=%f\n",sC33)); return kFALSE; } Double_t sC44=GetCovariance()[14]; if (sC44<=0 /*|| sC44>6e-5*/) { - if (n>fgkWARN) Warning("Invariant","fC44=%f\n",sC44); + if (n>fgkWARN) AliDebug(1,Form("fC44=%f\n",sC44)); return kFALSE; } @@ -346,7 +364,7 @@ Bool_t AliITStrackV2::Propagate(Double_t alp,Double_t xk) { if (!Invariant()) { Int_t n=GetNumberOfClusters(); - if (n>fgkWARN) AliWarning("Wrong invariant !"); + if (n>fgkWARN) AliDebug(1,"Wrong invariant !"); return kFALSE; } @@ -395,8 +413,6 @@ Bool_t AliITStrackV2::Improve(Double_t x0,Double_t xyz[3],Double_t ers[3]) { //------------------------------------------------------------------ //Store the initail track parameters - return kTRUE; //PH temporary switched off - Double_t x = GetX(); Double_t alpha = GetAlpha(); Double_t par[] = {GetY(),GetZ(),GetSnp(),GetTgl(),GetSigned1Pt()}; @@ -420,29 +436,31 @@ Bool_t AliITStrackV2::Improve(Double_t x0,Double_t xyz[3],Double_t ers[3]) { Double_t cs=TMath::Cos(GetAlpha()), sn=TMath::Sin(GetAlpha()); -//Double_t xv = xyz[0]*cs + xyz[1]*sn; // vertex + Double_t xv = xyz[0]*cs + xyz[1]*sn; // vertex Double_t yv =-xyz[0]*sn + xyz[1]*cs; // in the Double_t zv = xyz[2]; // local frame - Double_t dy = par[0] - yv, dz = par[1] - zv; - Double_t r2=GetX()*GetX() + dy*dy; + Double_t dx = x - xv, dy = par[0] - yv, dz = par[1] - zv; + Double_t r2=dx*dx + dy*dy; Double_t p2=(1.+ GetTgl()*GetTgl())/(GetSigned1Pt()*GetSigned1Pt()); Double_t beta2=p2/(p2 + GetMass()*GetMass()); x0*=TMath::Sqrt((1.+ GetTgl()*GetTgl())/(1.- GetSnp()*GetSnp())); Double_t theta2=14.1*14.1/(beta2*p2*1e6)*x0; //Double_t theta2=1.0259e-6*14*14/28/(beta2*p2)*x0*9.36*2.33; - Double_t cnv=GetBz()*kB2C; + Double_t bz=GetBz(); + Double_t cnv=bz*kB2C; + Double_t curv=GetC(bz); { - Double_t dummy = 4/r2 - GetC()*GetC(); + Double_t dummy = 4/r2 - curv*curv; if (dummy < 0) return kFALSE; - Double_t parp = 0.5*(GetC()*GetX() + dy*TMath::Sqrt(dummy)); + Double_t parp = 0.5*(curv*dx + dy*TMath::Sqrt(dummy)); Double_t sigma2p = theta2*(1.-GetSnp())*(1.+GetSnp())*(1. + GetTgl()*GetTgl()); Double_t ovSqr2 = 1./TMath::Sqrt(r2); Double_t tfact = ovSqr2*(1.-dy*ovSqr2)*(1.+dy*ovSqr2); sigma2p += cov[0]*tfact*tfact; sigma2p += ers[1]*ers[1]/r2; - sigma2p += 0.25*cov[14]*cnv*cnv*GetX()*GetX(); + sigma2p += 0.25*cov[14]*cnv*cnv*dx*dx; Double_t eps2p=sigma2p/(cov[5] + sigma2p); par[0] += cov[3]/(cov[5] + sigma2p)*(parp - GetSnp()); par[2] = eps2p*GetSnp() + (1 - eps2p)*parp; @@ -450,7 +468,7 @@ Bool_t AliITStrackV2::Improve(Double_t x0,Double_t xyz[3],Double_t ers[3]) { cov[3] *= eps2p; } { - Double_t parl=0.5*GetC()*dz/TMath::ASin(0.5*GetC()*TMath::Sqrt(r2)); + Double_t parl=0.5*curv*dz/TMath::ASin(0.5*curv*TMath::Sqrt(r2)); Double_t sigma2l=theta2; sigma2l += cov[2]/r2 + cov[0]*dy*dy*dz*dz/(r2*r2*r2); sigma2l += ers[2]*ers[2]/r2; @@ -470,7 +488,7 @@ Bool_t AliITStrackV2::Improve(Double_t x0,Double_t xyz[3],Double_t ers[3]) { return kTRUE; } -void AliITStrackV2::CookdEdx(Double_t low, Double_t up) { +void AliITStrackV2::CookdEdx(Double_t /*low*/, Double_t /*up*/) { //----------------------------------------------------------------- // This function calculates dE/dX within the "low" and "up" cuts. // Origin: Boris Batyunya, JINR, Boris.Batiounia@cern.ch @@ -504,15 +522,13 @@ void AliITStrackV2::CookdEdx(Double_t low, Double_t up) { } while (swap); - Double_t nl=low*nc, nu =up*nc; - Float_t sumamp = 0; - Float_t sumweight =0; + Double_t sumamp=0,sumweight=0; + Double_t weight[4]={1.,1.,0.,0.}; + if(nc==3) weight[1]=0.5; + else if(nc<3) weight[1]=0.; for (Int_t i=0; inu-1) weight = TMath::Max(nu-i,0.); - sumamp+= dedx[i]*weight; - sumweight+=weight; + sumamp+= dedx[i]*weight[i]; + sumweight+=weight[i]; } SetdEdx(sumamp/sumweight); } @@ -577,3 +593,178 @@ GetLocalXat(Double_t r,Double_t &xloc) const { return kTRUE; } + +//____________________________________________________________________________ +Bool_t AliITStrackV2::ImproveKalman(Double_t xyz[3],Double_t ers[3], const Double_t* xlMS, const Double_t* x2X0MS, Int_t nMS) +{ + // Substitute the state of the track (p_{k|k},C_{k|k}) at the k-th measumerent by its + // smoothed value from the k-th measurement + measurement at the vertex. + // Account for the MS on nMS layers at x-postions xlMS with x/x0 = x2X0MS + // p_{k|kv} = p_{k|k} + C_{k|k}*D^Tr_{k+1} B^{-1}_{k+1} ( vtx - D_{k+1}*p_{k|k}) + // C_{k|kv} = C_{k|k}*( I - D^Tr_{k+1} B^{-1}_{k+1} D_{k+1} C^Tr_{k|k}) + // + // where D_{k} = H_{k} F_{k} with H being the matrix converting the tracks parameters + // to measurements m_{k} = H_{k} p_{k} and F_{k} the matrix propagating the track between the + // the point k-1 and k: p_{k|k-1} = F_{k} p_{k-1|k-1} + // + // B_{k+1} = V_{k+1} + H_{k+1} C_{k+1|k} H^Tr_{k+1} with V_{k+1} being the error of the measurment + // at point k+1 (i.e. vertex), and C_{k+1|k} - error matrix extrapolated from k-th measurement to + // k+1 (vtx) and accounting for the MS inbetween + // + // H = {{1,0,0,0,0},{0,1,0,0,0}} + // + double covc[15], *cori = (double*) GetCovariance(),par[5] = {GetY(),GetZ(),GetSnp(),GetTgl(),GetSigned1Pt()}, + &c00=cori[0], + &c01=cori[1],&c11=cori[2], + &c02=cori[3],&c12=cori[4],&c22=cori[5], + &c03=cori[6],&c13=cori[7],&c23=cori[8],&c33=cori[9], + &c04=cori[10],&c14=cori[11],&c24=cori[12],&c34=cori[13],&c44=cori[14], + // for smoothed cov matrix + &cov00=covc[0], + &cov01=covc[1],&cov11=covc[2], + &cov02=covc[3],&cov12=covc[4],&cov22=covc[5], + &cov03=covc[6],&cov13=covc[7],&cov23=covc[8],&cov33=covc[9], + &cov04=covc[10],&cov14=covc[11],&cov24=covc[12],&cov34=covc[13],&cov44=covc[14]; + // + double x = GetX(), alpha = GetAlpha(); + // vertex in the track frame + double cs=TMath::Cos(alpha), sn=TMath::Sin(alpha); + double xv = xyz[0]*cs + xyz[1]*sn, yv =-xyz[0]*sn + xyz[1]*cs, zv = xyz[2]; + double dx = xv - GetX(); + if (TMath::Abs(dx)<=kAlmost0) return kTRUE; + // + double cnv=GetBz()*kB2C, x2r=cnv*par[4]*dx, f1=par[2], f2=f1+x2r; + if (TMath::Abs(f1) >= kAlmost1 || TMath::Abs(f2) >= kAlmost1) { + AliInfo(Form("Fail: %+e %+e",f1,f2)); + return kFALSE; + } + double r1=TMath::Sqrt((1.-f1)*(1.+f1)), r2=TMath::Sqrt((1.-f2)*(1.+f2)), dx2r=dx/(r1+r2); + // elements of matrix F_{k+1} (1s on diagonal) + double f02 = 2*dx2r, f04 = cnv*dx*dx2r, f13/*, f24 = cnv*dx*/; + if (TMath::Abs(x2r)<0.05) f13 = dx*r2+f2*(f1+f2)*dx2r; // see AliExternalTrackParam::PropagateTo + else { + double dy2dx = (f1+f2)/(r1+r2); + f13 = 2*TMath::ASin(0.5*TMath::Sqrt(1+dy2dx*dy2dx)*x2r)/(cnv*par[4]); + } + // elements of matrix D_{k+1} = H_{k+1} * F_{k+1} + // double d00 = 1., d11 = 1.; + double &d02 = f02, &d04 = f04, &d13 = f13; + // + // elements of matrix DC = D_{k+1}*C_{kk}^T + double dc00 = c00+c02*d02+c04*d04, dc10 = c01+c03*d13; + double dc01 = c01+c12*d02+c14*d04, dc11 = c11+c13*d13; + double dc02 = c02+c22*d02+c24*d04, dc12 = c12+c23*d13; + double dc03 = c03+c23*d02+c34*d04, dc13 = c13+c33*d13; + double dc04 = c04+c24*d02+c44*d04, dc14 = c14+c34*d13; + // + // difference between the vertex and the the track extrapolated to vertex + yv -= par[0] + par[2]*d02 + par[4]*d04; + zv -= par[1] + par[3]*d13; + // + // y,z part of the cov.matrix extrapolated to vtx (w/o MS contribution) + // C_{k+1,k} = H F_{k+1} C_{k,k} F^Tr_{k+1} H^Tr = D C D^Tr + double cv00 = dc00+dc02*d02+dc04*d04, cv01 = dc01+dc03*d13, cv11 = dc11+dc13*d13; + // + // add MS contribution layer by layer + double xCurr = x; + double p2Curr = par[2]; + // + // precalculated factors of MS contribution matrix: + double ms22t = (1. + par[3]*par[3]); + double ms33t = ms22t*ms22t; + double p34 = par[3]*par[4]; + double ms34t = p34*ms22t; + double ms44t = p34*p34; + // + double p2=(1.+ par[3]*par[3])/(par[4]*par[4]); + double beta2 = p2/(p2+GetMass()*GetMass()); + double theta2t = 14.1*14.1/(beta2*p2*1e6) * (1. + par[3]*par[3]); + // + // account for the MS in the layers between the last measurement and the vertex + for (int il=0;il= kAlmost1 || TMath::Abs(f2L) >= kAlmost1) { + AliInfo(Form("FailMS at step %d of %d: dfx:%e dxL:%e %e %e",il,nMS,dfx,dxL,f1L,f2L)); + return kFALSE; + } + double r1L=TMath::Sqrt((1.-f1L)*(1.+f1L)), r2L=TMath::Sqrt((1.-f2L)*(1.+f2L)), dx2rL=dxL/(r1L+r2L); + // elements of matrix for propagation from scatering layer to vertex + double f02L = 2*dx2rL, f04L = cnv*dxL*dx2rL, f13L/*, f24L = cnv*dxL*/; + if (TMath::Abs(x2rL)<0.05) f13L = dxL*r2L+f2L*(f1L+f2L)*dx2rL; // see AliExternalTrackParam::PropagateTo + else { + double dy2dxL = (f1L+f2L)/(r1L+r2L); + f13L = 2*TMath::ASin(0.5*TMath::Sqrt(1+dy2dxL*dy2dxL)*x2rL)/(cnv*par[4]); + } + // MS contribution matrix: + double theta2 = theta2t*TMath::Abs(x2X0MS[il]); + double ms22 = theta2*(1.-p2Curr)*(1.+p2Curr)*ms22t; + double ms33 = theta2*ms33t; + double ms34 = theta2*ms34t; + double ms44 = theta2*ms44t; + // + // add H F MS F^Tr H^Tr to cv + cv00 += f02L*f02L*ms22 + f04L*f04L*ms44; + cv01 += f04L*f13L*ms34; + cv11 += f13L*f13L*ms33; + } + // + // inverse of matrix B + double b11 = ers[1]*ers[1] + cv00; + double b00 = ers[2]*ers[2] + cv11; + double det = b11*b00 - cv01*cv01; + if (TMath::Abs(det)