]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCtrack.cxx
Setting the branch address to permit correct reading
[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 $Log$
18 Revision 1.19  2003/03/03 16:56:53  hristov
19 Corrections to obey coding conventions
20
21 Revision 1.18  2003/02/19 08:57:04  hristov
22 Control^M removed
23
24 Revision 1.17  2003/02/19 08:49:46  hristov
25 Track time measurement (S.Radomski)
26
27 Revision 1.16  2003/02/06 11:11:36  kowal2
28 Added a few get methods by Jiri Chudoba
29
30 Revision 1.15  2002/11/25 09:33:30  hristov
31 Tracking of secondaries (M.Ivanov)
32
33 Revision 1.14  2002/10/23 13:45:00  hristov
34 Fatal if no magnetic field set for the reconstruction (Y.Belikov)
35
36 Revision 1.13  2002/10/23 07:17:34  alibrary
37 Introducing Riostream.h
38
39 Revision 1.12  2002/10/14 14:57:43  hristov
40 Merging the VirtualMC branch to the main development branch (HEAD)
41
42 Revision 1.9.6.1  2002/10/11 08:34:48  hristov
43 Updating VirtualMC to v3-09-02
44
45 Revision 1.11  2002/07/19 07:34:42  kowal2
46 Logs added
47
48 */
49
50
51 //-----------------------------------------------------------------
52 //           Implementation of the TPC track class
53 //        This class is used by the AliTPCtracker class
54 //      Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
55 //-----------------------------------------------------------------
56
57 #include <Riostream.h>
58
59 #include "AliTPCtrack.h"
60 #include "AliCluster.h"
61 #include "AliBarrelTrack.h"
62
63 ClassImp(AliTPCtrack)
64
65 //_________________________________________________________________________
66 AliTPCtrack::AliTPCtrack(): AliKalmanTrack() 
67 {
68   fX = fP0 = fP1 = fP2 = fP3 = fP3 = fP4 = 0.0;
69   fAlpha = fdEdx = 0.0;
70   fNWrong = fNRotation = fNumber = 0;  // [SR, 01.04.2003]
71 }
72
73 //_________________________________________________________________________
74 AliTPCtrack::AliTPCtrack(UInt_t index, const Double_t xx[5],
75 const Double_t cc[15], Double_t xref, Double_t alpha) : AliKalmanTrack() {
76   //-----------------------------------------------------------------
77   // This is the main track constructor.
78   //-----------------------------------------------------------------
79   fX=xref;
80   fAlpha=alpha;
81   if (fAlpha<-TMath::Pi()) fAlpha += 2*TMath::Pi();
82   if (fAlpha>=TMath::Pi()) fAlpha -= 2*TMath::Pi();
83   fdEdx=0.;
84
85   fP0=xx[0]; fP1=xx[1]; fP2=xx[2]; fP3=xx[3]; fP4=xx[4];
86
87   fC00=cc[0];
88   fC10=cc[1];  fC11=cc[2];
89   fC20=cc[3];  fC21=cc[4];  fC22=cc[5];
90   fC30=cc[6];  fC31=cc[7];  fC32=cc[8];  fC33=cc[9];
91   fC40=cc[10]; fC41=cc[11]; fC42=cc[12]; fC43=cc[13]; fC44=cc[14];
92
93   fIndex[0]=index;
94   SetNumberOfClusters(1);
95
96 }
97
98 //_____________________________________________________________________________
99 AliTPCtrack::AliTPCtrack(const AliKalmanTrack& t,Double_t alpha) :
100 AliKalmanTrack(t) {
101   //-----------------------------------------------------------------
102   // Conversion AliKalmanTrack -> AliTPCtrack.
103   //-----------------------------------------------------------------
104   SetChi2(0.);
105   SetNumberOfClusters(0);
106
107   fdEdx  = 0.;
108   fAlpha = alpha;
109   if      (fAlpha < -TMath::Pi()) fAlpha += 2*TMath::Pi();
110   else if (fAlpha >= TMath::Pi()) fAlpha -= 2*TMath::Pi();
111
112   //Conversion of the track parameters
113   Double_t x,p[5]; t.GetExternalParameters(x,p);
114   fX=x;    x=GetConvConst();
115   fP0=p[0]; 
116   fP1=p[1]; 
117   fP3=p[3];
118   fP4=p[4]/x; 
119   fP2=fP4*fX - p[2];
120
121   //Conversion of the covariance matrix
122   Double_t c[15]; t.GetExternalCovariance(c);
123   c[10]/=x; c[11]/=x; c[12]/=x; c[13]/=x; c[14]/=x*x;
124
125   Double_t c22=fX*fX*c[14] - 2*fX*c[12] + c[5];
126   Double_t c32=fX*c[13] - c[8];
127   Double_t c20=fX*c[10] - c[3], c21=fX*c[11] - c[4], c42=fX*c[14] - c[12];
128   
129   fC00=c[0 ];
130   fC10=c[1 ];   fC11=c[2 ];
131   fC20=c20;     fC21=c21;     fC22=c22;
132   fC30=c[6 ];   fC31=c[7 ];   fC32=c32;   fC33=c[9 ];
133   fC40=c[10];   fC41=c[11];   fC42=c42;   fC43=c[13]; fC44=c[14];
134
135 }
136
137 //_____________________________________________________________________________
138 AliTPCtrack::AliTPCtrack(const AliTPCtrack& t) : AliKalmanTrack(t) {
139   //-----------------------------------------------------------------
140   // This is a track copy constructor.
141   //-----------------------------------------------------------------
142   fX=t.fX;
143   fAlpha=t.fAlpha;
144   fdEdx=t.fdEdx;
145
146   fP0=t.fP0; fP1=t.fP1; fP2=t.fP2; fP3=t.fP3; fP4=t.fP4;
147
148   fC00=t.fC00;
149   fC10=t.fC10;  fC11=t.fC11;
150   fC20=t.fC20;  fC21=t.fC21;  fC22=t.fC22;
151   fC30=t.fC30;  fC31=t.fC31;  fC32=t.fC32;  fC33=t.fC33;
152   fC40=t.fC40;  fC41=t.fC41;  fC42=t.fC42;  fC43=t.fC43;  fC44=t.fC44;
153
154   Int_t n=GetNumberOfClusters();
155   for (Int_t i=0; i<n; i++) fIndex[i]=t.fIndex[i];
156 }
157 //_____________________________________________________________________________
158
159 void  AliTPCtrack::GetBarrelTrack(AliBarrelTrack *track) {
160   //
161   // Create a Barrel Track out of this track
162   // Current track is propagated to the reference plane
163   // by the tracker
164   //
165   // [SR, 01.04.2003]
166   
167   if (!track) return;
168   Double_t xr, vec[5], cov[15];
169
170   track->SetLabel(GetLabel());
171   track->SetX(fX, fAlpha);
172   track->SetNClusters(GetNumberOfClusters(), GetChi2());
173   track->SetTime(fIntegratedTime, fIntegratedLength);
174
175   track->SetMass(fMass);
176   track->SetdEdX(fdEdx);
177
178   track->SetNWrongClusters(fNWrong);
179   track->SetNRotate(fNRotation);
180
181   GetExternalParameters(xr, vec);
182   track->SetStateVector(vec);
183   
184   GetExternalCovariance(cov);
185   track->SetCovarianceMatrix(cov);
186
187 }
188 //_____________________________________________________________________________
189 Int_t AliTPCtrack::Compare(const TObject *o) const {
190   //-----------------------------------------------------------------
191   // This function compares tracks according to the their curvature
192   //-----------------------------------------------------------------
193   AliTPCtrack *t=(AliTPCtrack*)o;
194   //Double_t co=TMath::Abs(t->Get1Pt());
195   //Double_t c =TMath::Abs(Get1Pt());
196   Double_t co=t->GetSigmaY2()*t->GetSigmaZ2();
197   Double_t c =GetSigmaY2()*GetSigmaZ2();
198   if (c>co) return 1;
199   else if (c<co) return -1;
200   return 0;
201 }
202
203 //_____________________________________________________________________________
204 void AliTPCtrack::GetExternalCovariance(Double_t cc[15]) const {
205   //-------------------------------------------------------------------------
206   // This function returns an external representation of the covriance matrix.
207   //   (See comments in AliTPCtrack.h about external track representation)
208   //-------------------------------------------------------------------------
209   Double_t a=GetConvConst();
210
211   Double_t c22=fX*fX*fC44-2*fX*fC42+fC22;
212   Double_t c32=fX*fC43-fC32;
213   Double_t c20=fX*fC40-fC20, c21=fX*fC41-fC21, c42=fX*fC44-fC42;
214   
215   cc[0 ]=fC00;
216   cc[1 ]=fC10;   cc[2 ]=fC11;
217   cc[3 ]=c20;    cc[4 ]=c21;    cc[5 ]=c22;
218   cc[6 ]=fC30;   cc[7 ]=fC31;   cc[8 ]=c32;   cc[9 ]=fC33; 
219   cc[10]=fC40*a; cc[11]=fC41*a; cc[12]=c42*a; cc[13]=fC43*a; cc[14]=fC44*a*a;
220
221 }
222
223 //_____________________________________________________________________________
224 Double_t AliTPCtrack::GetPredictedChi2(const AliCluster *c) const 
225 {
226   //-----------------------------------------------------------------
227   // This function calculates a predicted chi2 increment.
228   //-----------------------------------------------------------------
229   Double_t r00=c->GetSigmaY2(), r01=0., r11=c->GetSigmaZ2();
230   r00+=fC00; r01+=fC10; r11+=fC11;
231
232   Double_t det=r00*r11 - r01*r01;
233   if (TMath::Abs(det) < 1.e-10) {
234     Int_t n=GetNumberOfClusters();
235     if (n>4) cerr<<n<<" AliKalmanTrack warning: Singular matrix !\n";
236     return 1e10;
237   }
238   Double_t tmp=r00; r00=r11; r11=tmp; r01=-r01;
239   
240   Double_t dy=c->GetY() - fP0, dz=c->GetZ() - fP1;
241   
242   return (dy*r00*dy + 2*r01*dy*dz + dz*r11*dz)/det;
243 }
244
245 //_____________________________________________________________________________
246 Int_t AliTPCtrack::PropagateTo(Double_t xk,Double_t x0,Double_t rho) {
247   //-----------------------------------------------------------------
248   // This function propagates a track to a reference plane x=xk.
249   //-----------------------------------------------------------------
250   if (TMath::Abs(fP4*xk - fP2) >= 0.9) {
251     //    Int_t n=GetNumberOfClusters();
252     //if (n>4) cerr<<n<<" AliTPCtrack warning: Propagation failed !\n";
253     return 0;
254   }
255   
256   // old position for time [SR, GSI 17.02.2003]
257   Double_t oldX = fX;
258   Double_t oldY = fP0;
259   Double_t oldZ = fP1;
260   //
261
262   Double_t x1=fX, x2=x1+(xk-x1), dx=x2-x1, y1=fP0, z1=fP1;
263   Double_t c1=fP4*x1 - fP2, r1=sqrt(1.- c1*c1);
264   Double_t c2=fP4*x2 - fP2, r2=sqrt(1.- c2*c2);
265   
266   fP0 += dx*(c1+c2)/(r1+r2);
267   fP1 += dx*(c1+c2)/(c1*r2 + c2*r1)*fP3;
268
269   //f = F - 1
270   Double_t rr=r1+r2, cc=c1+c2, xx=x1+x2;
271   Double_t f02=-dx*(2*rr + cc*(c1/r1 + c2/r2))/(rr*rr);
272   Double_t f04= dx*(rr*xx + cc*(c1*x1/r1+c2*x2/r2))/(rr*rr);
273   Double_t cr=c1*r2+c2*r1;
274   Double_t f12=-dx*fP3*(2*cr + cc*(c2*c1/r1-r1 + c1*c2/r2-r2))/(cr*cr);
275   Double_t f13= dx*cc/cr; 
276   Double_t f14=dx*fP3*(cr*xx-cc*(r1*x2-c2*c1*x1/r1+r2*x1-c1*c2*x2/r2))/(cr*cr);
277
278   //b = C*ft
279   Double_t b00=f02*fC20 + f04*fC40, b01=f12*fC20 + f14*fC40 + f13*fC30;
280   Double_t b10=f02*fC21 + f04*fC41, b11=f12*fC21 + f14*fC41 + f13*fC31;
281   Double_t b20=f02*fC22 + f04*fC42, b21=f12*fC22 + f14*fC42 + f13*fC32;
282   Double_t b30=f02*fC32 + f04*fC43, b31=f12*fC32 + f14*fC43 + f13*fC33;
283   Double_t b40=f02*fC42 + f04*fC44, b41=f12*fC42 + f14*fC44 + f13*fC43;
284   
285   //a = f*b = f*C*ft
286   Double_t a00=f02*b20+f04*b40,a01=f02*b21+f04*b41,a11=f12*b21+f14*b41+f13*b31;
287
288   //F*C*Ft = C + (a + b + bt)
289   fC00 += a00 + 2*b00;
290   fC10 += a01 + b01 + b10; 
291   fC20 += b20;
292   fC30 += b30;
293   fC40 += b40;
294   fC11 += a11 + 2*b11;
295   fC21 += b21; 
296   fC31 += b31; 
297   fC41 += b41; 
298
299   fX=x2;
300
301   //Multiple scattering******************
302   Double_t d=sqrt((x1-fX)*(x1-fX)+(y1-fP0)*(y1-fP0)+(z1-fP1)*(z1-fP1));
303   Double_t p2=(1.+ GetTgl()*GetTgl())/(Get1Pt()*Get1Pt());
304   Double_t beta2=p2/(p2 + GetMass()*GetMass());
305   //Double_t theta2=14.1*14.1/(beta2*p2*1e6)*d/x0*rho;
306   Double_t theta2=1.0259e-6*10*10/20/(beta2*p2)*d*rho;
307
308   Double_t ey=fP4*fX - fP2, ez=fP3;
309   Double_t xz=fP4*ez, zz1=ez*ez+1, xy=fP2+ey;
310
311   fC22 += (2*ey*ez*ez*fP2+1-ey*ey+ez*ez+fP2*fP2*ez*ez)*theta2;
312   fC32 += ez*zz1*xy*theta2;
313   fC33 += zz1*zz1*theta2;
314   fC42 += xz*ez*xy*theta2;
315   fC43 += xz*zz1*theta2;
316   fC44 += xz*xz*theta2;
317
318   //Energy losses************************
319   Double_t dE=0.153e-3/beta2*(log(5940*beta2/(1-beta2)) - beta2)*d*rho;
320   if (x1 < x2) dE=-dE;
321   cc=fP4;
322   fP4*=(1.- sqrt(p2+GetMass()*GetMass())/p2*dE);
323   fP2+=fX*(fP4-cc);
324
325   // Integrated Time [SR, GSI, 17.02.2003]
326   if (IsStartedTimeIntegral()) {
327     Double_t l2 = (fX-oldX)*(fX-oldX)+(fP0-oldY)*(fP0-oldY)+(fP1-oldZ)*(fP1-oldZ);
328     AddTimeStep(TMath::Sqrt(l2));
329   }
330   //
331
332   return 1;
333 }
334
335 //_____________________________________________________________________________
336 Int_t AliTPCtrack::PropagateToVertex(Double_t x0,Double_t rho) 
337 {
338   //-----------------------------------------------------------------
339   // This function propagates tracks to the "vertex".
340   //-----------------------------------------------------------------
341   Double_t c=fP4*fX - fP2;
342   Double_t tgf=-fP2/(fP4*fP0 + sqrt(1-c*c));
343   Double_t snf=tgf/sqrt(1.+ tgf*tgf);
344   Double_t xv=(fP2+snf)/fP4;
345   return PropagateTo(xv,x0,rho);
346 }
347
348 //_____________________________________________________________________________
349 Int_t AliTPCtrack::Update(const AliCluster *c, Double_t chisq, UInt_t index) {
350   //-----------------------------------------------------------------
351   // This function associates a cluster with this track.
352   //-----------------------------------------------------------------
353
354   // update the number of wrong SR[20.03.2003]
355   Int_t absLabel = TMath::Abs(GetLabel());
356   if ( (c->GetLabel(0) != absLabel) && 
357        (c->GetLabel(0) != absLabel) &&
358        (c->GetLabel(0) != absLabel)) fNWrong++;
359   //
360
361   Double_t r00=c->GetSigmaY2(), r01=0., r11=c->GetSigmaZ2();
362   r00+=fC00; r01+=fC10; r11+=fC11;
363   Double_t det=r00*r11 - r01*r01;
364   Double_t tmp=r00; r00=r11/det; r11=tmp/det; r01=-r01/det;
365
366   Double_t k00=fC00*r00+fC10*r01, k01=fC00*r01+fC10*r11;
367   Double_t k10=fC10*r00+fC11*r01, k11=fC10*r01+fC11*r11;
368   Double_t k20=fC20*r00+fC21*r01, k21=fC20*r01+fC21*r11;
369   Double_t k30=fC30*r00+fC31*r01, k31=fC30*r01+fC31*r11;
370   Double_t k40=fC40*r00+fC41*r01, k41=fC40*r01+fC41*r11;
371
372   Double_t dy=c->GetY() - fP0, dz=c->GetZ() - fP1;
373   Double_t cur=fP4 + k40*dy + k41*dz, eta=fP2 + k20*dy + k21*dz;
374   if (TMath::Abs(cur*fX-eta) >= 0.9) {
375     //    Int_t n=GetNumberOfClusters();
376     //if (n>4) cerr<<n<<" AliTPCtrack warning: Filtering failed !\n";
377     return 0;
378   }
379
380   fP0 += k00*dy + k01*dz;
381   fP1 += k10*dy + k11*dz;
382   fP2  = eta;
383   fP3 += k30*dy + k31*dz;
384   fP4  = cur;
385
386   Double_t c01=fC10, c02=fC20, c03=fC30, c04=fC40;
387   Double_t c12=fC21, c13=fC31, c14=fC41;
388
389   fC00-=k00*fC00+k01*fC10; fC10-=k00*c01+k01*fC11;
390   fC20-=k00*c02+k01*c12;   fC30-=k00*c03+k01*c13;
391   fC40-=k00*c04+k01*c14; 
392
393   fC11-=k10*c01+k11*fC11;
394   fC21-=k10*c02+k11*c12;   fC31-=k10*c03+k11*c13;
395   fC41-=k10*c04+k11*c14; 
396
397   fC22-=k20*c02+k21*c12;   fC32-=k20*c03+k21*c13;
398   fC42-=k20*c04+k21*c14; 
399
400   fC33-=k30*c03+k31*c13;
401   fC43-=k40*c03+k41*c13; 
402
403   fC44-=k40*c04+k41*c14; 
404
405   Int_t n=GetNumberOfClusters();
406   fIndex[n]=index;
407   SetNumberOfClusters(n+1);
408   SetChi2(GetChi2()+chisq);
409
410   return 1;
411 }
412
413 //_____________________________________________________________________________
414 Int_t AliTPCtrack::Rotate(Double_t alpha)
415 {
416   //-----------------------------------------------------------------
417   // This function rotates this track.
418   //-----------------------------------------------------------------
419
420   if (alpha != 0) fNRotation++;  // [SR, 01.04.2003]
421
422   fAlpha += alpha;
423   if (fAlpha<-TMath::Pi()) fAlpha += 2*TMath::Pi();
424   if (fAlpha>=TMath::Pi()) fAlpha -= 2*TMath::Pi();
425   
426   Double_t x1=fX, y1=fP0;
427   Double_t ca=cos(alpha), sa=sin(alpha);
428   Double_t r1=fP4*fX - fP2;
429   
430   fX = x1*ca + y1*sa;
431   fP0=-x1*sa + y1*ca;
432   fP2=fP2*ca + (fP4*y1 + sqrt(1.- r1*r1))*sa;
433   
434   Double_t r2=fP4*fX - fP2;
435   if (TMath::Abs(r2) >= 0.99999) {
436     Int_t n=GetNumberOfClusters();
437     if (n>4) cerr<<n<<" AliTPCtrack warning: Rotation failed !\n";
438     return 0;
439   }
440   
441   Double_t y0=fP0 + sqrt(1.- r2*r2)/fP4;
442   if ((fP0-y0)*fP4 >= 0.) {
443     Int_t n=GetNumberOfClusters();
444     if (n>4) cerr<<n<<" AliTPCtrack warning: Rotation failed !!!\n";
445     return 0;
446   }
447
448   //f = F - 1
449   Double_t f00=ca-1,    f24=(y1 - r1*x1/sqrt(1.- r1*r1))*sa, 
450            f20=fP4*sa,  f22=(ca + sa*r1/sqrt(1.- r1*r1))-1;
451
452   //b = C*ft
453   Double_t b00=fC00*f00, b02=fC00*f20+fC40*f24+fC20*f22;
454   Double_t b10=fC10*f00, b12=fC10*f20+fC41*f24+fC21*f22;
455   Double_t b20=fC20*f00, b22=fC20*f20+fC42*f24+fC22*f22;
456   Double_t b30=fC30*f00, b32=fC30*f20+fC43*f24+fC32*f22;
457   Double_t b40=fC40*f00, b42=fC40*f20+fC44*f24+fC42*f22;
458
459   //a = f*b = f*C*ft
460   Double_t a00=f00*b00, a02=f00*b02, a22=f20*b02+f24*b42+f22*b22;
461
462   // *** Double_t dy2=fCyy;
463
464   //F*C*Ft = C + (a + b + bt)
465   fC00 += a00 + 2*b00;
466   fC10 += b10;
467   fC20 += a02+b20+b02;
468   fC30 += b30;
469   fC40 += b40;
470   fC21 += b12;
471   fC32 += b32;
472   fC22 += a22 + 2*b22;
473   fC42 += b42; 
474
475   // *** fCyy+=dy2*sa*sa*r1*r1/(1.- r1*r1);
476   // *** fCzz+=d2y*sa*sa*fT*fT/(1.- r1*r1);
477
478   return 1;
479 }
480
481 void AliTPCtrack::ResetCovariance() {
482   //------------------------------------------------------------------
483   //This function makes a track forget its history :)  
484   //------------------------------------------------------------------
485
486   fC00*=10.;
487   fC10=0.;  fC11*=10.;
488   fC20=0.;  fC21=0.;  fC22*=10.;
489   fC30=0.;  fC31=0.;  fC32=0.;  fC33*=10.;
490   fC40=0.;  fC41=0.;  fC42=0.;  fC43=0.;  fC44*=10.;
491
492 }
493
494 ////////////////////////////////////////////////////////////////////////
495 Double_t AliTPCtrack::Phi() const {
496 //
497 //
498 //
499   Double_t phi =  TMath::ASin(GetSnp()) + fAlpha;
500   if (phi<0) phi+=2*TMath::Pi();
501   if (phi>=2*TMath::Pi()) phi-=2*TMath::Pi();
502   return phi;
503 }
504 ////////////////////////////////////////////////////////////////////////
505