]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCtrack.cxx
New design of tracking classes (Yu.Belikov)
[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 //-----------------------------------------------------------------
17 //           Implementation of the TPC track class
18 //
19 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
20 //-----------------------------------------------------------------
21
22 #include <iostream.h>
23
24 #include "AliCluster.h"
25 #include "AliTPCtrack.h"
26 #include "AliTPCcluster.h"
27 #include "AliTPCClustersRow.h"
28 #include "AliTPCClustersArray.h"
29
30 ClassImp(AliTPCtrack)
31 //_________________________________________________________________________
32 AliTPCtrack::AliTPCtrack(UInt_t index, const Double_t xx[5],
33 const Double_t cc[15], Double_t xref, Double_t alpha) : AliKalmanTrack() {
34   //-----------------------------------------------------------------
35   // This is the main track constructor.
36   //-----------------------------------------------------------------
37   fX=xref;
38   fAlpha=alpha;
39   if (fAlpha<-TMath::Pi()) fAlpha += 2*TMath::Pi();
40   if (fAlpha>=TMath::Pi()) fAlpha -= 2*TMath::Pi();
41   fdEdx=0.;
42
43   fP0=xx[0]; fP1=xx[1]; fP2=xx[2]; fP3=xx[3]; fP4=xx[4];
44
45   fC00=cc[0];
46   fC10=cc[1];  fC11=cc[2];
47   fC20=cc[3];  fC21=cc[4];  fC22=cc[5];
48   fC30=cc[6];  fC31=cc[7];  fC32=cc[8];  fC33=cc[9];
49   fC40=cc[10]; fC41=cc[11]; fC42=cc[12]; fC43=cc[13]; fC44=cc[14];
50
51   fIndex[0]=index;
52   SetNumberOfClusters(1);
53 }
54
55 //_____________________________________________________________________________
56 AliTPCtrack::AliTPCtrack(const AliTPCtrack& t) : AliKalmanTrack(t) {
57   //-----------------------------------------------------------------
58   // This is a track copy constructor.
59   //-----------------------------------------------------------------
60   fX=t.fX;
61   fAlpha=t.fAlpha;
62   fdEdx=t.fdEdx;
63
64   fP0=t.fP0; fP1=t.fP1; fP2=t.fP2; fP3=t.fP3; fP4=t.fP4;
65
66   fC00=t.fC00;
67   fC10=t.fC10;  fC11=t.fC11;
68   fC20=t.fC20;  fC21=t.fC21;  fC22=t.fC22;
69   fC30=t.fC30;  fC31=t.fC31;  fC32=t.fC32;  fC33=t.fC33;
70   fC40=t.fC40;  fC41=t.fC41;  fC42=t.fC42;  fC43=t.fC43;  fC44=t.fC44;
71
72   Int_t n=GetNumberOfClusters();
73   for (Int_t i=0; i<n; i++) fIndex[i]=t.fIndex[i];
74 }
75
76 //_____________________________________________________________________________
77 Int_t AliTPCtrack::Compare(const TObject *o) const {
78   //-----------------------------------------------------------------
79   // This function compares tracks according to the their curvature
80   //-----------------------------------------------------------------
81   AliTPCtrack *t=(AliTPCtrack*)o;
82   Double_t co=TMath::Abs(t->Get1Pt());
83   Double_t c =TMath::Abs(Get1Pt());
84   if (c>co) return 1;
85   else if (c<co) return -1;
86   return 0;
87 }
88
89 //_____________________________________________________________________________
90 void AliTPCtrack::GetExternalCovariance(Double_t cc[15]) const {
91   //-------------------------------------------------------------------------
92   // This function returns an external representation of the covriance matrix.
93   //   (See comments in AliTPCtrack.h about external track representation)
94   //-------------------------------------------------------------------------
95   Double_t a=kConversionConstant;
96
97   Double_t c22=fX*fX*fC33-2*fX*fC32+fC22;
98   Double_t c42=fX*fC43-fC42;
99   Double_t c20=fX*fC30-fC20, c21=fX*fC31-fC21, c32=fX*fC33-fC32;
100   
101   cc[0 ]=fC00;
102   cc[1 ]=fC10;   cc[2 ]=fC11;
103   cc[3 ]=c20;    cc[4 ]=c21;    cc[5 ]=c22;
104   cc[6 ]=fC40;   cc[7 ]=fC41;   cc[8 ]=c42;   cc[9 ]=fC44; 
105   cc[10]=fC30*a; cc[11]=fC31*a; cc[12]=c32*a; cc[13]=fC43*a; cc[14]=fC33*a*a;
106
107 }
108
109 //_____________________________________________________________________________
110 Double_t AliTPCtrack::GetPredictedChi2(const AliCluster *c) const 
111 {
112   //-----------------------------------------------------------------
113   // This function calculates a predicted chi2 increment.
114   //-----------------------------------------------------------------
115   Double_t r00=c->GetSigmaY2(), r01=0., r11=c->GetSigmaZ2();
116   r00+=fC00; r01+=fC10; r11+=fC11;
117
118   Double_t det=r00*r11 - r01*r01;
119   if (TMath::Abs(det) < 1.e-10) {
120     Int_t n=GetNumberOfClusters();
121     if (n>4) cerr<<n<<" AliKalmanTrack warning: Singular matrix !\n";
122     return 1e10;
123   }
124   Double_t tmp=r00; r00=r11; r11=tmp; r01=-r01;
125   
126   Double_t dy=c->GetY() - fP0, dz=c->GetZ() - fP1;
127   
128   return (dy*r00*dy + 2*r01*dy*dz + dz*r11*dz)/det;
129 }
130
131 //_____________________________________________________________________________
132 Int_t AliTPCtrack::PropagateTo(Double_t xk,Double_t x0,Double_t rho,Double_t pm)
133 {
134   //-----------------------------------------------------------------
135   // This function propagates a track to a reference plane x=xk.
136   //-----------------------------------------------------------------
137   if (TMath::Abs(fP3*xk - fP2) >= 0.99999) {
138     Int_t n=GetNumberOfClusters();
139     if (n>4) cerr<<n<<" AliTPCtrack warning: Propagation failed !\n";
140     return 0;
141   }
142
143   Double_t x1=fX, x2=x1+(xk-x1), dx=x2-x1, y1=fP0, z1=fP1;
144   Double_t c1=fP3*x1 - fP2, r1=sqrt(1.- c1*c1);
145   Double_t c2=fP3*x2 - fP2, r2=sqrt(1.- c2*c2);
146   
147   fP0 += dx*(c1+c2)/(r1+r2);
148   fP1 += dx*(c1+c2)/(c1*r2 + c2*r1)*fP4;
149
150   //f = F - 1
151   Double_t rr=r1+r2, cc=c1+c2, xx=x1+x2;
152   Double_t f02=-dx*(2*rr + cc*(c1/r1 + c2/r2))/(rr*rr);
153   Double_t f03= dx*(rr*xx + cc*(c1*x1/r1+c2*x2/r2))/(rr*rr);
154   Double_t cr=c1*r2+c2*r1;
155   Double_t f12=-dx*fP4*(2*cr + cc*(c2*c1/r1-r1 + c1*c2/r2-r2))/(cr*cr);
156   Double_t f13=dx*fP4*(cr*xx-cc*(r1*x2-c2*c1*x1/r1+r2*x1-c1*c2*x2/r2))/(cr*cr);
157   Double_t f14= dx*cc/cr; 
158
159   //b = C*ft
160   Double_t b00=f02*fC20 + f03*fC30, b01=f12*fC20 + f13*fC30 + f14*fC40;
161   Double_t b10=f02*fC21 + f03*fC31, b11=f12*fC21 + f13*fC31 + f14*fC41;
162   Double_t b20=f02*fC22 + f03*fC32, b21=f12*fC22 + f13*fC32 + f14*fC42;
163   Double_t b30=f02*fC32 + f03*fC33, b31=f12*fC32 + f13*fC33 + f14*fC43;
164   Double_t b40=f02*fC42 + f03*fC43, b41=f12*fC42 + f13*fC43 + f14*fC44;
165   
166   //a = f*b = f*C*ft
167   Double_t a00=f02*b20+f03*b30,a01=f02*b21+f03*b31,a11=f12*b21+f13*b31+f14*b41;
168
169   //F*C*Ft = C + (a + b + bt)
170   fC00 += a00 + 2*b00;
171   fC10 += a01 + b01 + b10; 
172   fC20 += b20;
173   fC30 += b30;
174   fC40 += b40;
175   fC11 += a11 + 2*b11;
176   fC21 += b21; 
177   fC31 += b31; 
178   fC41 += b41; 
179
180   fX=x2;
181
182   //Multiple scattering******************
183   Double_t d=sqrt((x1-fX)*(x1-fX)+(y1-fP0)*(y1-fP0)+(z1-fP1)*(z1-fP1));
184   Double_t p2=(1.+ GetTgl()*GetTgl())/(Get1Pt()*Get1Pt());
185   Double_t beta2=p2/(p2 + pm*pm);
186   Double_t theta2=14.1*14.1/(beta2*p2*1e6)*d/x0*rho;
187   //Double_t theta2=1.0259e-6*10*10/20/(beta2*p2)*d*rho;
188
189   Double_t ey=fP3*fX - fP2, ez=fP4;
190   Double_t xz=fP3*ez, zz1=ez*ez+1, xy=fP2+ey;
191
192   fC33 += xz*xz*theta2;
193   fC32 += xz*ez*xy*theta2;
194   fC43 += xz*zz1*theta2;
195   fC22 += (2*ey*ez*ez*fP2+1-ey*ey+ez*ez+fP2*fP2*ez*ez)*theta2;
196   fC42 += ez*zz1*xy*theta2;
197   fC44 += zz1*zz1*theta2;
198
199   //Energy losses************************
200   Double_t dE=0.153e-3/beta2*(log(5940*beta2/(1-beta2)) - beta2)*d*rho;
201   if (x1 < x2) dE=-dE;
202   cc=fP3;
203   fP3*=(1.- sqrt(p2+pm*pm)/p2*dE);
204   fP2+=fX*(fP3-cc);
205
206   return 1;
207 }
208
209 //_____________________________________________________________________________
210 Int_t AliTPCtrack::PropagateToVertex(Double_t x0,Double_t rho,Double_t pm) 
211 {
212   //-----------------------------------------------------------------
213   // This function propagates tracks to the "vertex".
214   //-----------------------------------------------------------------
215   Double_t c=fP3*fX - fP2;
216   Double_t tgf=-fP2/(fP3*fP0 + sqrt(1-c*c));
217   Double_t snf=tgf/sqrt(1.+ tgf*tgf);
218   Double_t xv=(fP2+snf)/fP3;
219   return PropagateTo(xv,x0,rho,pm);
220 }
221
222 //_____________________________________________________________________________
223 Int_t AliTPCtrack::Update(const AliCluster *c, Double_t chisq, UInt_t index) {
224   //-----------------------------------------------------------------
225   // This function associates a cluster with this track.
226   //-----------------------------------------------------------------
227   Double_t r00=c->GetSigmaY2(), r01=0., r11=c->GetSigmaZ2();
228   r00+=fC00; r01+=fC10; r11+=fC11;
229   Double_t det=r00*r11 - r01*r01;
230   Double_t tmp=r00; r00=r11/det; r11=tmp/det; r01=-r01/det;
231
232   Double_t k00=fC00*r00+fC10*r01, k01=fC00*r01+fC10*r11;
233   Double_t k10=fC10*r00+fC11*r01, k11=fC10*r01+fC11*r11;
234   Double_t k20=fC20*r00+fC21*r01, k21=fC20*r01+fC21*r11;
235   Double_t k30=fC30*r00+fC31*r01, k31=fC30*r01+fC31*r11;
236   Double_t k40=fC40*r00+fC41*r01, k41=fC40*r01+fC41*r11;
237
238   Double_t dy=c->GetY() - fP0, dz=c->GetZ() - fP1;
239   Double_t cur=fP3 + k30*dy + k31*dz, eta=fP2 + k20*dy + k21*dz;
240   if (TMath::Abs(cur*fX-eta) >= 0.99999) {
241     Int_t n=GetNumberOfClusters();
242     if (n>4) cerr<<n<<" AliTPCtrack warning: Filtering failed !\n";
243     return 0;
244   }
245
246   fP0 += k00*dy + k01*dz;
247   fP1 += k10*dy + k11*dz;
248   fP2  = eta;
249   fP3  = cur;
250   fP4 += k40*dy + k41*dz;
251
252   Double_t c01=fC10, c02=fC20, c03=fC30, c04=fC40;
253   Double_t c12=fC21, c13=fC31, c14=fC41;
254
255   fC00-=k00*fC00+k01*fC10; fC10-=k00*c01+k01*fC11;
256   fC20-=k00*c02+k01*c12;   fC30-=k00*c03+k01*c13;
257   fC40-=k00*c04+k01*c14; 
258
259   fC11-=k10*c01+k11*fC11;
260   fC21-=k10*c02+k11*c12;   fC31-=k10*c03+k11*c13;
261   fC41-=k10*c04+k11*c14; 
262
263   fC22-=k20*c02+k21*c12;   fC32-=k20*c03+k21*c13;
264   fC42-=k20*c04+k21*c14; 
265
266   fC33-=k30*c03+k31*c13;
267   fC43-=k30*c04+k31*c14; 
268
269   fC44-=k40*c04+k41*c14; 
270
271   Int_t n=GetNumberOfClusters();
272   fIndex[n]=index;
273   SetNumberOfClusters(n+1);
274   SetChi2(GetChi2()+chisq);
275
276   return 1;
277 }
278
279 //_____________________________________________________________________________
280 Int_t AliTPCtrack::Rotate(Double_t alpha)
281 {
282   //-----------------------------------------------------------------
283   // This function rotates this track.
284   //-----------------------------------------------------------------
285   fAlpha += alpha;
286   if (fAlpha<-TMath::Pi()) fAlpha += 2*TMath::Pi();
287   if (fAlpha>=TMath::Pi()) fAlpha -= 2*TMath::Pi();
288   
289   Double_t x1=fX, y1=fP0;
290   Double_t ca=cos(alpha), sa=sin(alpha);
291   Double_t r1=fP3*fX - fP2;
292   
293   fX = x1*ca + y1*sa;
294   fP0=-x1*sa + y1*ca;
295   fP2=fP2*ca + (fP3*y1 + sqrt(1.- r1*r1))*sa;
296   
297   Double_t r2=fP3*fX - fP2;
298   if (TMath::Abs(r2) >= 0.99999) {
299     Int_t n=GetNumberOfClusters();
300     if (n>4) cerr<<n<<" AliTPCtrack warning: Rotation failed !\n";
301     return 0;
302   }
303   
304   Double_t y0=fP0 + sqrt(1.- r2*r2)/fP3;
305   if ((fP0-y0)*fP3 >= 0.) {
306     Int_t n=GetNumberOfClusters();
307     if (n>4) cerr<<n<<" AliTPCtrack warning: Rotation failed !!!\n";
308     return 0;
309   }
310
311   //f = F - 1
312   Double_t f00=ca-1,    f23=(y1 - r1*x1/sqrt(1.- r1*r1))*sa, 
313            f20=fP3*sa,  f22=(ca + sa*r1/sqrt(1.- r1*r1))-1;
314
315   //b = C*ft
316   Double_t b00=fC00*f00, b02=fC00*f20+fC30*f23+fC20*f22;
317   Double_t b10=fC10*f00, b12=fC10*f20+fC31*f23+fC21*f22;
318   Double_t b20=fC20*f00, b22=fC20*f20+fC32*f23+fC22*f22;
319   Double_t b30=fC30*f00, b32=fC30*f20+fC33*f23+fC32*f22;
320   Double_t b40=fC40*f00, b42=fC40*f20+fC43*f23+fC42*f22;
321
322   //a = f*b = f*C*ft
323   Double_t a00=f00*b00, a02=f00*b02, a22=f20*b02+f23*b32+f22*b22;
324
325   // *** Double_t dy2=fCyy;
326
327   //F*C*Ft = C + (a + b + bt)
328   fC00 += a00 + 2*b00;
329   fC10 += b10;
330   fC20 += a02+b20+b02;
331   fC30 += b30;
332   fC40 += b40;
333   fC21 += b12;
334   fC32 += b32;
335   fC22 += a22 + 2*b22;
336   fC42 += b42; 
337
338   // *** fCyy+=dy2*sa*sa*r1*r1/(1.- r1*r1);
339   // *** fCzz+=d2y*sa*sa*fT*fT/(1.- r1*r1);
340
341   return 1;
342 }
343
344 //_____________________________________________________________________________
345 void AliTPCtrack::CookLabel(AliTPCClustersArray *ca) {
346   //-----------------------------------------------------------------
347   // This function cooks the track label. If label<0, this track is fake.
348   //-----------------------------------------------------------------
349   Int_t n=GetNumberOfClusters();
350   Int_t *lb=new Int_t[n];
351   Int_t *mx=new Int_t[n];
352   AliTPCcluster **clusters=new AliTPCcluster*[n];
353
354   Int_t i;
355   Int_t sec,row,ncl;
356   for (i=0; i<n; i++) {
357      lb[i]=mx[i]=0;
358      GetCluster(i,sec,row,ncl);
359      AliTPCClustersRow *clrow=ca->GetRow(sec,row);
360      clusters[i]=(AliTPCcluster*)(*clrow)[ncl];      
361   }
362   
363   Int_t lab=123456789;
364   for (i=0; i<n; i++) {
365     AliTPCcluster *c=clusters[i];
366     lab=TMath::Abs(c->GetLabel(0));
367     Int_t j;
368     for (j=0; j<n; j++)
369       if (lb[j]==lab || mx[j]==0) break;
370     lb[j]=lab;
371     (mx[j])++;
372   }
373   
374   Int_t max=0;
375   for (i=0; i<n; i++) 
376     if (mx[i]>max) {max=mx[i]; lab=lb[i];}
377     
378   for (i=0; i<n; i++) {
379     AliTPCcluster *c=clusters[i];
380     if (TMath::Abs(c->GetLabel(1)) == lab ||
381         TMath::Abs(c->GetLabel(2)) == lab ) max++;
382   }
383   
384   SetLabel(-lab);
385   if (1.-Float_t(max)/n <= 0.10) {
386     //Int_t tail=Int_t(0.08*fN);
387      Int_t tail=14;
388      max=0;
389      for (i=1; i<=tail; i++) {
390        AliTPCcluster *c=clusters[n-i];
391        if (lab == TMath::Abs(c->GetLabel(0)) ||
392            lab == TMath::Abs(c->GetLabel(1)) ||
393            lab == TMath::Abs(c->GetLabel(2))) max++;
394      }
395      if (max >= Int_t(0.5*tail)) SetLabel(lab);
396   }
397 /*
398 if (lab==111)
399    for (i=0; i<n; i++) {
400      AliTPCcluster *c=clusters[i];
401      cerr<<i<<' '<<c->GetLabel(0)<<endl;
402    }
403 */
404   delete[] lb;
405   delete[] mx;
406   delete[] clusters;
407 }