]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCtrack.cxx
Corrected FMD data in the ESD (C.Holm)
[u/mrichter/AliRoot.git] / TPC / AliTPCtrack.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17
18 //-----------------------------------------------------------------
19 //           Implementation of the TPC track class
20 //        This class is used by the AliTPCtracker class
21 //      Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
22 //-----------------------------------------------------------------
23
24 #include <Riostream.h>
25
26 #include "AliTPCtrack.h"
27 #include "AliCluster.h"
28 #include "AliESDtrack.h"
29
30 ClassImp(AliTPCtrack)
31
32 //_________________________________________________________________________
33 AliTPCtrack::AliTPCtrack(): AliKalmanTrack() 
34 {
35   //-------------------------------------------------
36   // default constructor
37   //-------------------------------------------------
38   fX = fP0 = fP1 = fP2 = fP3 = fP3 = fP4 = 0.0;
39   fAlpha = fdEdx = 0.0;
40   fNumber = 0;  // [SR, 01.04.2003]
41   for (Int_t i=0; i<3;i++) fKinkIndexes[i]=0;
42 }
43
44 //_________________________________________________________________________
45
46
47
48 AliTPCtrack::AliTPCtrack(UInt_t index, const Double_t xx[5],
49 const Double_t cc[15], Double_t xref, Double_t alpha) : AliKalmanTrack() {
50   //-----------------------------------------------------------------
51   // This is the main track constructor.
52   //-----------------------------------------------------------------
53   fX=xref;
54   fAlpha=alpha;
55   if (fAlpha<-TMath::Pi()) fAlpha += 2*TMath::Pi();
56   if (fAlpha>=TMath::Pi()) fAlpha -= 2*TMath::Pi();
57   fdEdx=0.;
58
59   fP0=xx[0]; fP1=xx[1]; fP2=xx[2]; fP3=xx[3]; fP4=xx[4];
60
61   SaveLocalConvConst();
62
63   fC00=cc[0];
64   fC10=cc[1];  fC11=cc[2];
65   fC20=cc[3];  fC21=cc[4];  fC22=cc[5];
66   fC30=cc[6];  fC31=cc[7];  fC32=cc[8];  fC33=cc[9];
67   fC40=cc[10]; fC41=cc[11]; fC42=cc[12]; fC43=cc[13]; fC44=cc[14];
68
69   fIndex[0]=index;
70   SetNumberOfClusters(1);
71   //
72   //MI
73   fSdEdx      = 0;
74   fNFoundable = 0;
75   fBConstrain = 0;
76   fLastPoint  = 0;
77   fFirstPoint = 0;
78   fRemoval    = 0;
79   fTrackType  = 0;
80   fLab2       = 0;
81   for (Int_t i=0; i<3;i++) fKinkIndexes[i]=0;
82   for (Int_t i=0; i<3;i++) fV0Indexes[i]=0;
83 }
84
85 //_____________________________________________________________________________
86 AliTPCtrack::AliTPCtrack(const AliESDtrack& t) : AliKalmanTrack() {
87   //-----------------------------------------------------------------
88   // Conversion AliESDtrack -> AliTPCtrack.
89   //-----------------------------------------------------------------
90   SetNumberOfClusters(t.GetTPCclusters(fIndex));
91   SetLabel(t.GetLabel());
92   SetMass(t.GetMass());
93   for (Int_t i=0; i<3;i++) fKinkIndexes[i]=t.GetKinkIndex(i);
94   for (Int_t i=0; i<3;i++) fV0Indexes[i]=t.GetV0Index(i);
95
96   fdEdx  = t.GetTPCsignal();
97   fAlpha = t.GetAlpha();
98   if      (fAlpha < -TMath::Pi()) fAlpha += 2*TMath::Pi();
99   else if (fAlpha >= TMath::Pi()) fAlpha -= 2*TMath::Pi();
100
101   //Conversion of the track parameters
102   Double_t x,p[5]; t.GetExternalParameters(x,p);
103   Double_t c[15]; t.GetExternalCovariance(c);
104
105   fX=x;    
106   fP0=p[0];
107   fP1=p[1];    SaveLocalConvConst();
108   fP3=p[3];    x=GetLocalConvConst();
109   fP4=p[4]/x;
110   fP2=fP4*fX - p[2];
111
112   //Conversion of the covariance matrix
113   c[10]/=x; c[11]/=x; c[12]/=x; c[13]/=x; c[14]/=x*x;
114
115   Double_t c22=fX*fX*c[14] - 2*fX*c[12] + c[5];
116   Double_t c32=fX*c[13] - c[8];
117   Double_t c20=fX*c[10] - c[3], c21=fX*c[11] - c[4], c42=fX*c[14] - c[12];
118
119   fC00=c[0 ];
120   fC10=c[1 ];   fC11=c[2 ];
121   fC20=c20;     fC21=c21;     fC22=c22;
122   fC30=c[6 ];   fC31=c[7 ];   fC32=c32;   fC33=c[9 ];
123   fC40=c[10];   fC41=c[11];   fC42=c42;   fC43=c[13]; fC44=c[14];
124
125   if ((t.GetStatus()&AliESDtrack::kTIME) == 0) return;
126   StartTimeIntegral();
127   Double_t times[10]; t.GetIntegratedTimes(times); SetIntegratedTimes(times);
128   SetIntegratedLength(t.GetIntegratedLength());
129   //
130   //MI
131   fSdEdx      = 0;
132   fNFoundable = 0;
133   fBConstrain = 0;
134   fLastPoint  = 0;
135   fFirstPoint = 0;
136   fRemoval    = 0;
137   fTrackType  = 0;
138   fLab2       = 0;
139   //  SetFakeRatio(t.GetTPCFakeRatio());
140 }
141
142 //_____________________________________________________________________________
143 AliTPCtrack::AliTPCtrack(const AliTPCtrack& t) : AliKalmanTrack(t) {
144   //-----------------------------------------------------------------
145   // This is a track copy constructor.
146   //-----------------------------------------------------------------
147   fX=t.fX;
148   fAlpha=t.fAlpha;
149   fdEdx=t.fdEdx;
150
151   fP0=t.fP0; fP1=t.fP1; fP2=t.fP2; fP3=t.fP3; fP4=t.fP4;
152
153   fC00=t.fC00;
154   fC10=t.fC10;  fC11=t.fC11;
155   fC20=t.fC20;  fC21=t.fC21;  fC22=t.fC22;
156   fC30=t.fC30;  fC31=t.fC31;  fC32=t.fC32;  fC33=t.fC33;
157   fC40=t.fC40;  fC41=t.fC41;  fC42=t.fC42;  fC43=t.fC43;  fC44=t.fC44;
158
159   //Int_t n=GetNumberOfClusters();
160   for (Int_t i=0; i<kMaxRow; i++) fIndex[i]=t.fIndex[i];
161   //
162   //MI 
163   fSdEdx      = t.fSdEdx;
164   fNFoundable = t.fNFoundable;
165   fBConstrain = t.fBConstrain;
166   fLastPoint  = t.fLastPoint;
167   fFirstPoint = t.fFirstPoint;
168   fRemoval    = t.fRemoval ;
169   fTrackType  = t.fTrackType;
170   fLab2       = t.fLab2;
171   for (Int_t i=0; i<3;i++) fKinkIndexes[i]=t.fKinkIndexes[i];
172   for (Int_t i=0; i<3;i++) fV0Indexes[i]=t.fV0Indexes[i];
173 }
174
175 //_____________________________________________________________________________
176 Int_t AliTPCtrack::Compare(const TObject *o) const {
177   //-----------------------------------------------------------------
178   // This function compares tracks according to the their curvature
179   //-----------------------------------------------------------------
180   AliTPCtrack *t=(AliTPCtrack*)o;
181   //Double_t co=TMath::Abs(t->Get1Pt());
182   //Double_t c =TMath::Abs(Get1Pt());
183   Double_t co=t->GetSigmaY2()*t->GetSigmaZ2();
184   Double_t c =GetSigmaY2()*GetSigmaZ2();
185   if (c>co) return 1;
186   else if (c<co) return -1;
187   return 0;
188 }
189
190 //_____________________________________________________________________________
191 void AliTPCtrack::GetExternalCovariance(Double_t cc[15]) const {
192  //-------------------------------------------------------------------------
193   // This function returns an external representation of the covriance matrix.
194   //   (See comments in AliTPCtrack.h about external track representation)
195   //-------------------------------------------------------------------------
196   Double_t a=GetLocalConvConst();
197
198   Double_t c22=fX*fX*fC44-2*fX*fC42+fC22;
199   Double_t c32=fX*fC43-fC32;
200   Double_t c20=fX*fC40-fC20, c21=fX*fC41-fC21, c42=fX*fC44-fC42;
201   
202   cc[0 ]=fC00;
203   cc[1 ]=fC10;   cc[2 ]=fC11;
204   cc[3 ]=c20;    cc[4 ]=c21;    cc[5 ]=c22;
205   cc[6 ]=fC30;   cc[7 ]=fC31;   cc[8 ]=c32;   cc[9 ]=fC33; 
206   cc[10]=fC40*a; cc[11]=fC41*a; cc[12]=c42*a; cc[13]=fC43*a; cc[14]=fC44*a*a;
207
208 }
209
210 //_____________________________________________________________________________
211 Double_t AliTPCtrack::GetPredictedChi2(const AliCluster *c) const 
212 {
213   //-----------------------------------------------------------------
214   // This function calculates a predicted chi2 increment.
215   //-----------------------------------------------------------------
216   Double_t r00=c->GetSigmaY2(), r01=0., r11=c->GetSigmaZ2();
217   r00+=fC00; r01+=fC10; r11+=fC11;
218
219   Double_t det=r00*r11 - r01*r01;
220   if (TMath::Abs(det) < 1.e-10) {
221     Int_t n=GetNumberOfClusters();
222     if (n>4) cerr<<n<<" AliKalmanTrack warning: Singular matrix !\n";
223     return 1e10;
224   }
225   Double_t tmp=r00; r00=r11; r11=tmp; r01=-r01;
226   
227   Double_t dy=c->GetY() - fP0, dz=c->GetZ() - fP1;
228   
229   return (dy*r00*dy + 2*r01*dy*dz + dz*r11*dz)/det;
230 }
231
232 Double_t AliTPCtrack::GetYat(Double_t xk) const {
233 //-----------------------------------------------------------------
234 // This function calculates the Y-coordinate of a track at the plane x=xk.
235 //-----------------------------------------------------------------
236   if (TMath::Abs(fP4*fX - fP2)>0.9999) return 1e10; //patch 01 jan 06
237     Double_t c1=fP4*fX - fP2, r1=TMath::Sqrt(1.- c1*c1);
238     Double_t c2=fP4*xk - fP2;
239     if (c2*c2>0.99999) {
240        Int_t n=GetNumberOfClusters();
241        if (n>4) cerr<<n<<"AliTPCtrack::GetYat: can't evaluate the y-coord !\n";
242        return 1e10;
243     } 
244     Double_t r2=TMath::Sqrt(1.- c2*c2);
245     return fP0 + (xk-fX)*(c1+c2)/(r1+r2);
246 }
247
248 //_____________________________________________________________________________
249 Int_t AliTPCtrack::PropagateTo(Double_t xk,Double_t /*x0*/,Double_t rho) {
250   //-----------------------------------------------------------------
251   // This function propagates a track to a reference plane x=xk.
252   //-----------------------------------------------------------------
253   if (TMath::Abs(fP4*xk - fP2) >= 0.9) {
254     //    Int_t n=GetNumberOfClusters();
255     //if (n>4) cerr<<n<<" AliTPCtrack warning: Propagation failed !\n";
256     return 0;
257   }
258   Double_t lcc=GetLocalConvConst();  
259
260   // old position for time [SR, GSI 17.02.2003]
261   Double_t oldX = fX;
262   Double_t oldY = fP0;
263   Double_t oldZ = fP1;
264   //
265
266   Double_t x1=fX, x2=x1+(xk-x1), dx=x2-x1, y1=fP0, z1=fP1;
267   Double_t c1=fP4*x1 - fP2, r1=sqrt(1.- c1*c1);
268   Double_t c2=fP4*x2 - fP2, r2=sqrt(1.- c2*c2);
269   
270   fP0 += dx*(c1+c2)/(r1+r2);
271   fP1 += dx*(c1+c2)/(c1*r2 + c2*r1)*fP3;
272
273   //f = F - 1
274   Double_t rr=r1+r2, cc=c1+c2, xx=x1+x2;
275   Double_t f02=-dx*(2*rr + cc*(c1/r1 + c2/r2))/(rr*rr);
276   Double_t f04= dx*(rr*xx + cc*(c1*x1/r1+c2*x2/r2))/(rr*rr);
277   Double_t cr=c1*r2+c2*r1;
278   Double_t f12=-dx*fP3*(2*cr + cc*(c2*c1/r1-r1 + c1*c2/r2-r2))/(cr*cr);
279   Double_t f13= dx*cc/cr; 
280   Double_t f14=dx*fP3*(cr*xx-cc*(r1*x2-c2*c1*x1/r1+r2*x1-c1*c2*x2/r2))/(cr*cr);
281
282   //b = C*ft
283   Double_t b00=f02*fC20 + f04*fC40, b01=f12*fC20 + f14*fC40 + f13*fC30;
284   Double_t b10=f02*fC21 + f04*fC41, b11=f12*fC21 + f14*fC41 + f13*fC31;
285   Double_t b20=f02*fC22 + f04*fC42, b21=f12*fC22 + f14*fC42 + f13*fC32;
286   Double_t b30=f02*fC32 + f04*fC43, b31=f12*fC32 + f14*fC43 + f13*fC33;
287   Double_t b40=f02*fC42 + f04*fC44, b41=f12*fC42 + f14*fC44 + f13*fC43;
288   
289   //a = f*b = f*C*ft
290   Double_t a00=f02*b20+f04*b40,a01=f02*b21+f04*b41,a11=f12*b21+f14*b41+f13*b31;
291
292   //F*C*Ft = C + (a + b + bt)
293   fC00 += a00 + 2*b00;
294   fC10 += a01 + b01 + b10; 
295   fC20 += b20;
296   fC30 += b30;
297   fC40 += b40;
298   fC11 += a11 + 2*b11;
299   fC21 += b21; 
300   fC31 += b31; 
301   fC41 += b41; 
302
303   fX=x2;
304
305   //Change of the magnetic field *************
306   SaveLocalConvConst();
307   cc=fP4;
308   fP4*=lcc/GetLocalConvConst();
309   fP2+=fX*(fP4-cc);
310
311   //Multiple scattering ******************
312   Double_t d=sqrt((x1-fX)*(x1-fX)+(y1-fP0)*(y1-fP0)+(z1-fP1)*(z1-fP1));
313   Double_t p2=(1.+ GetTgl()*GetTgl())/(Get1Pt()*Get1Pt());
314   Double_t beta2=p2/(p2 + GetMass()*GetMass());
315   beta2 = TMath::Min(beta2,0.99999999999);
316   //Double_t theta2=14.1*14.1/(beta2*p2*1e6)*d/x0*rho;
317   Double_t theta2=1.0259e-6*10*10/20/(beta2*p2)*d*rho;
318
319   Double_t ey=fP4*fX - fP2, ez=fP3;
320   Double_t xz=fP4*ez, zz1=ez*ez+1, xy=fP2+ey;
321     
322   fC22 += (2*ey*ez*ez*fP2+1-ey*ey+ez*ez+fP2*fP2*ez*ez)*theta2;
323   fC32 += ez*zz1*xy*theta2;
324   fC33 += zz1*zz1*theta2;
325   fC42 += xz*ez*xy*theta2;
326   fC43 += xz*zz1*theta2;
327   fC44 += xz*xz*theta2;
328   /*
329   //
330   //MI coeficints
331   Double_t dc22 = (1-ey*ey+xz*xz*fX*fX)*theta2;
332   Double_t dc32 = (xz*fX*zz1)*theta2;
333   Double_t dc33 = (zz1*zz1)*theta2;
334   Double_t dc42 = (xz*fX*xz)*theta2;
335   Double_t dc43 = (zz1*xz)*theta2;
336   Double_t dc44 = (xz*xz)*theta2; 
337   fC22 += dc22;
338   fC32 += dc32;
339   fC33 += dc33;
340   fC42 += dc42;
341   fC43 += dc43;
342   fC44 += dc44;
343   */
344   //Energy losses ************************
345   Double_t dE=0.153e-3/beta2*(log(5940*beta2/(1-beta2)) - beta2)*d*rho;
346   if (x1 < x2) dE=-dE;
347   cc=fP4;
348
349   //Double_t E = sqrt(p2+GetMass()*GetMass());
350   //Double_t mifac  = TMath::Sqrt(1.+dE*dE/p2+2*E*dE/p2)-1;
351   //Double_t belfac = E*dE/p2;
352                                //
353   fP4*=(1.- sqrt(p2+GetMass()*GetMass())/p2*dE);
354   fP2+=fX*(fP4-cc);
355
356   // Integrated Time [SR, GSI, 17.02.2003]
357  if (x1 < x2)
358  if (IsStartedTimeIntegral()) {
359     Double_t l2 = (fX-oldX)*(fX-oldX)+(fP0-oldY)*(fP0-oldY)+(fP1-oldZ)*(fP1-oldZ);
360     AddTimeStep(TMath::Sqrt(l2));
361   }
362   //
363
364   return 1;
365 }
366
367 //_____________________________________________________________________________
368 Int_t AliTPCtrack::PropagateToVertex(Double_t x0,Double_t rho) 
369 {
370   //-----------------------------------------------------------------
371   // This function propagates tracks to the "vertex".
372   //-----------------------------------------------------------------
373   Double_t c=fP4*fX - fP2;
374   Double_t tgf=-fP2/(fP4*fP0 + sqrt(1-c*c));
375   Double_t snf=tgf/sqrt(1.+ tgf*tgf);
376   Double_t xv=(fP2+snf)/fP4;
377   return PropagateTo(xv,x0,rho);
378 }
379
380 //_____________________________________________________________________________
381 Int_t AliTPCtrack::Update(const AliCluster *c, Double_t chisq, UInt_t index) {
382   //-----------------------------------------------------------------
383   // This function associates a cluster with this track.
384   //-----------------------------------------------------------------
385   Double_t r00=c->GetSigmaY2(), r01=0., r11=c->GetSigmaZ2();
386   r00+=fC00; r01+=fC10; r11+=fC11;
387   Double_t det=r00*r11 - r01*r01;
388   Double_t tmp=r00; r00=r11/det; r11=tmp/det; r01=-r01/det;
389
390   Double_t k00=fC00*r00+fC10*r01, k01=fC00*r01+fC10*r11;
391   Double_t k10=fC10*r00+fC11*r01, k11=fC10*r01+fC11*r11;
392   Double_t k20=fC20*r00+fC21*r01, k21=fC20*r01+fC21*r11;
393   Double_t k30=fC30*r00+fC31*r01, k31=fC30*r01+fC31*r11;
394   Double_t k40=fC40*r00+fC41*r01, k41=fC40*r01+fC41*r11;
395
396   Double_t dy=c->GetY() - fP0, dz=c->GetZ() - fP1;
397   Double_t cur=fP4 + k40*dy + k41*dz, eta=fP2 + k20*dy + k21*dz;
398   if (TMath::Abs(cur*fX-eta) >= 0.9) {
399     //    Int_t n=GetNumberOfClusters();
400     //if (n>4) cerr<<n<<" AliTPCtrack warning: Filtering failed !\n";
401     return 0;
402   }
403
404   fP0 += k00*dy + k01*dz;
405   fP1 += k10*dy + k11*dz;
406   fP2  = eta;
407   fP3 += k30*dy + k31*dz;
408   fP4  = cur;
409
410   Double_t c01=fC10, c02=fC20, c03=fC30, c04=fC40;
411   Double_t c12=fC21, c13=fC31, c14=fC41;
412
413   fC00-=k00*fC00+k01*fC10; fC10-=k00*c01+k01*fC11;
414   fC20-=k00*c02+k01*c12;   fC30-=k00*c03+k01*c13;
415   fC40-=k00*c04+k01*c14; 
416
417   fC11-=k10*c01+k11*fC11;
418   fC21-=k10*c02+k11*c12;   fC31-=k10*c03+k11*c13;
419   fC41-=k10*c04+k11*c14; 
420
421   fC22-=k20*c02+k21*c12;   fC32-=k20*c03+k21*c13;
422   fC42-=k20*c04+k21*c14; 
423
424   fC33-=k30*c03+k31*c13;
425   fC43-=k40*c03+k41*c13; 
426
427   fC44-=k40*c04+k41*c14; 
428
429   Int_t n=GetNumberOfClusters();
430   fIndex[n]=index;
431   SetNumberOfClusters(n+1);
432   SetChi2(GetChi2()+chisq);
433
434   return 1;
435 }
436
437 //_____________________________________________________________________________
438 Int_t AliTPCtrack::Rotate(Double_t alpha)
439 {
440   //-----------------------------------------------------------------
441   // This function rotates this track.
442   //-----------------------------------------------------------------
443   fAlpha += alpha;
444   if (fAlpha<-TMath::Pi()) fAlpha += 2*TMath::Pi();
445   if (fAlpha>=TMath::Pi()) fAlpha -= 2*TMath::Pi();
446   
447   Double_t x1=fX, y1=fP0;
448   Double_t ca=cos(alpha), sa=sin(alpha);
449   Double_t r1=fP4*fX - fP2;
450   
451   if (TMath::Abs(r1)>=0.9999) return 0; //patch 01 jan 06
452
453   fX = x1*ca + y1*sa;
454   fP0=-x1*sa + y1*ca;
455   fP2=fP2*ca + (fP4*y1 + sqrt(1.- r1*r1))*sa;
456   
457   Double_t r2=fP4*fX - fP2;
458   if (TMath::Abs(r2) >= 0.99999) {
459     Int_t n=GetNumberOfClusters();
460     if (n>4) cerr<<n<<" AliTPCtrack warning: Rotation failed !\n";
461     return 0;
462   }
463   
464   Double_t y0=fP0 + sqrt(1.- r2*r2)/fP4;
465   if ((fP0-y0)*fP4 >= 0.) {
466     Int_t n=GetNumberOfClusters();
467     if (n>4) cerr<<n<<" AliTPCtrack warning: Rotation failed !!!\n";
468     return 0;
469   }
470
471   //f = F - 1
472   Double_t f00=ca-1,    f24=(y1 - r1*x1/sqrt(1.- r1*r1))*sa, 
473            f20=fP4*sa,  f22=(ca + sa*r1/sqrt(1.- r1*r1))-1;
474
475   //b = C*ft
476   Double_t b00=fC00*f00, b02=fC00*f20+fC40*f24+fC20*f22;
477   Double_t b10=fC10*f00, b12=fC10*f20+fC41*f24+fC21*f22;
478   Double_t b20=fC20*f00, b22=fC20*f20+fC42*f24+fC22*f22;
479   Double_t b30=fC30*f00, b32=fC30*f20+fC43*f24+fC32*f22;
480   Double_t b40=fC40*f00, b42=fC40*f20+fC44*f24+fC42*f22;
481
482   //a = f*b = f*C*ft
483   Double_t a00=f00*b00, a02=f00*b02, a22=f20*b02+f24*b42+f22*b22;
484
485   // *** Double_t dy2=fCyy;
486
487   //F*C*Ft = C + (a + b + bt)
488   fC00 += a00 + 2*b00;
489   fC10 += b10;
490   fC20 += a02+b20+b02;
491   fC30 += b30;
492   fC40 += b40;
493   fC21 += b12;
494   fC32 += b32;
495   fC22 += a22 + 2*b22;
496   fC42 += b42; 
497
498   // *** fCyy+=dy2*sa*sa*r1*r1/(1.- r1*r1);
499   // *** fCzz+=d2y*sa*sa*fT*fT/(1.- r1*r1);
500
501   return 1;
502 }
503
504 void AliTPCtrack::ResetCovariance() {
505   //------------------------------------------------------------------
506   //This function makes a track forget its history :)  
507   //------------------------------------------------------------------
508
509   fC00*=10.;
510   fC10=0.;  fC11*=10.;
511   fC20=0.;  fC21=0.;  fC22*=10.;
512   fC30=0.;  fC31=0.;  fC32=0.;  fC33*=10.;
513   fC40=0.;  fC41=0.;  fC42=0.;  fC43=0.;  fC44*=10.;
514
515 }
516
517 ////////////////////////////////////////////////////////////////////////
518 Double_t AliTPCtrack::Phi() const {
519 //
520 //
521 //
522   Double_t phi =  TMath::ASin(GetSnp()) + fAlpha;
523   if (phi<0) phi+=2*TMath::Pi();
524   if (phi>=2*TMath::Pi()) phi-=2*TMath::Pi();
525   return phi;
526 }
527 ////////////////////////////////////////////////////////////////////////
528
529
530
531 ////////////////////////////////////////////////////////////////////////
532 // MI ADDITION
533
534 Float_t AliTPCtrack::Density(Int_t row0, Int_t row1)
535 {
536   //
537   // calculate cluster density
538   Int_t good  = 0;
539   Int_t found = 0;
540   //if (row0<fFirstPoint) row0 = fFirstPoint;
541   if (row1>fLastPoint) row1 = fLastPoint;
542
543   
544   for (Int_t i=row0;i<=row1;i++){ 
545     //    Int_t index = fClusterIndex[i];
546     Int_t index = fIndex[i];
547     if (index!=-1)  good++;
548     if (index>0)    found++;
549   }
550   Float_t density=0;
551   if (good>0) density = Float_t(found)/Float_t(good);
552   return density;
553 }
554
555
556 Float_t AliTPCtrack::Density2(Int_t row0, Int_t row1)
557 {
558   //
559   // calculate cluster density
560   Int_t good  = 0;
561   Int_t found = 0;
562   //  
563   for (Int_t i=row0;i<=row1;i++){     
564     Int_t index = fIndex[i];
565     if (index!=-1)  good++;
566     if (index>0)    found++;
567   }
568   Float_t density=0;
569   if (good>0) density = Float_t(found)/Float_t(good);
570   return density;
571 }
572
573
574 Double_t AliTPCtrack::GetZat0() const
575 {
576   //
577   // return virtual z - supposing that x = 0
578   if (TMath::Abs(fP2)>1) return 0;
579   if (TMath::Abs(fX*fP4-fP2)>1) return 0;
580   Double_t vz = fP1+fP3/fP4*(asin(-fP2)-asin(fX*fP4-fP2));
581   return vz;
582 }
583
584
585 Double_t AliTPCtrack::GetD(Double_t x, Double_t y) const {
586   //------------------------------------------------------------------
587   // This function calculates the transverse impact parameter
588   // with respect to a point with global coordinates (x,y)
589   //------------------------------------------------------------------
590   //Double_t xt=fX, yt=fP0;
591
592   Double_t sn=TMath::Sin(fAlpha), cs=TMath::Cos(fAlpha);
593   Double_t a = x*cs + y*sn;
594   y = -x*sn + y*cs; x=a;
595   //
596   Double_t r  = TMath::Abs(1/fP4);
597   Double_t x0 = TMath::Abs(fP2*r);
598   Double_t y0 = fP0;
599   y0= fP0+TMath::Sqrt(1-(fP4*fX-fP2)*(fP4*fX-fP2))/fP4;
600   
601   Double_t  delta = TMath::Sqrt((x-x0)*(x-x0)+(y-y0)*(y-y0));  
602   //  Double_t  delta = TMath::Sqrt(TMath::Abs(x*x-2*x0*x+x0*x0+ y*y-2*y*y0+y0*y0));
603   delta -= TMath::Abs(r);
604   return delta;  
605 }
606
607 //
608 //
609
610 void  AliTPCtrack::UpdatePoints()
611 {
612   //--------------------------------------------------
613   //calculates first ,amx dens and last points
614   //--------------------------------------------------
615   Float_t density[160];
616   for (Int_t i=0;i<160;i++) density[i]=-1.;
617   fPoints[0]= 160;
618   fPoints[1] = -1;
619   //
620   Int_t ngood=0;
621   Int_t undeff=0;
622   Int_t nall =0;
623   Int_t range=20;
624   for (Int_t i=0;i<160;i++){
625     Int_t last = i-range;
626     if (nall<range) nall++;
627     if (last>=0){
628       if (fIndex[last]>0&& (fIndex[last]&0x8000)==0) ngood--;
629       if (fIndex[last]==-1) undeff--;
630     }
631     if (fIndex[i]>0&& (fIndex[i]&0x8000)==0)   ngood++;
632     if (fIndex[i]==-1) undeff++;
633     if (nall==range &&undeff<range/2) density[i-range/2] = Float_t(ngood)/Float_t(nall-undeff);
634   }
635   Float_t maxdens=0;
636   Int_t indexmax =0;
637   for (Int_t i=0;i<160;i++){
638     if (density[i]<0) continue;
639     if (density[i]>maxdens){
640       maxdens=density[i];
641       indexmax=i;
642     }
643   }
644   //
645   //max dens point
646   fPoints[3] = maxdens;
647   fPoints[1] = indexmax;
648   //
649   // last point
650   for (Int_t i=indexmax;i<160;i++){
651     if (density[i]<0) continue;
652     if (density[i]<maxdens/2.) {
653       break;
654     }
655     fPoints[2]=i;
656   }
657   //
658   // first point
659   for (Int_t i=indexmax;i>0;i--){
660     if (density[i]<0) continue;
661     if (density[i]<maxdens/2.) {
662       break;
663     }
664     fPoints[0]=i;
665   }
666   //
667 }