]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TOF/AliTOFtrack.cxx
AliAODParticle -> AliVAODParticle; AliAODStdParticle -> AliAODParticle
[u/mrichter/AliRoot.git] / TOF / AliTOFtrack.cxx
CommitLineData
74ea065c 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#include <Riostream.h>
19#include <TObject.h>
20#include "AliTOFGeometry.h"
21#include "AliTOFtrack.h"
22#include "AliESDtrack.h"
23
24ClassImp(AliTOFtrack)
25
26//_____________________________________________________________________________
27AliTOFtrack::AliTOFtrack(const AliTOFtrack& t) : AliKalmanTrack(t) {
28 //
29 // Copy constructor.
30 //
31
32 SetSeedIndex(t.GetSeedIndex());
33 SetLabel(t.GetLabel());
34 fSeedLab=t.GetSeedLabel();
35 SetChi2(t.GetChi2());
36
37 fAlpha=t.fAlpha;
38 fX=t.fX;
39
40 fY=t.fY; fZ=t.fZ; fE=t.fE; fT=t.fT; fC=t.fC;
41
42 fCyy=t.fCyy;
43 fCzy=t.fCzy; fCzz=t.fCzz;
44 fCey=t.fCey; fCez=t.fCez; fCee=t.fCee;
45 fCty=t.fCty; fCtz=t.fCtz; fCte=t.fCte; fCtt=t.fCtt;
46 fCcy=t.fCcy; fCcz=t.fCcz; fCce=t.fCce; fCct=t.fCct; fCcc=t.fCcc;
47
48
49}
50
51//_____________________________________________________________________________
52AliTOFtrack::AliTOFtrack(const AliESDtrack& t)
53 :AliKalmanTrack() {
54 //
55 // Constructor from AliESDtrack
56 //
57
58 SetSeedIndex(-1);
59 SetLabel(t.GetLabel());
60 SetChi2(0.);
61 SetMass(t.GetMass());
62
63 fAlpha = t.GetAlpha();
64 if (fAlpha < -TMath::Pi()) fAlpha += 2*TMath::Pi();
65 else if (fAlpha >= TMath::Pi()) fAlpha -= 2*TMath::Pi();
66 Double_t x, p[5]; t.GetExternalParameters(x,p);
67
68 fX=x;
69
70 x = GetConvConst();
71
72 fY=p[0];
73 fZ=p[1];
74 fT=p[3];
75 fC=p[4]/x;
76 fE=fC*fX - p[2];
77
78 //Conversion of the covariance matrix
79 Double_t c[15]; t.GetExternalCovariance(c);
80
81 c[10]/=x; c[11]/=x; c[12]/=x; c[13]/=x; c[14]/=x*x;
82
83 Double_t c22=fX*fX*c[14] - 2*fX*c[12] + c[5];
84 Double_t c32=fX*c[13] - c[8];
85 Double_t c20=fX*c[10] - c[3], c21=fX*c[11] - c[4], c42=fX*c[14] - c[12];
86
87 fCyy=c[0 ];
88 fCzy=c[1 ]; fCzz=c[2 ];
89 fCey=c20; fCez=c21; fCee=c22;
90 fCty=c[6 ]; fCtz=c[7 ]; fCte=c32; fCtt=c[9 ];
91 fCcy=c[10]; fCcz=c[11]; fCce=c42; fCct=c[13]; fCcc=c[14];
92
93 if ((t.GetStatus()&AliESDtrack::kTIME) == 0) return;
94 StartTimeIntegral();
95 Double_t times[10]; t.GetIntegratedTimes(times); SetIntegratedTimes(times);
96 SetIntegratedLength(t.GetIntegratedLength());
97
98
99}
100//____________________________________________________________________________
101void AliTOFtrack::GetExternalParameters(Double_t& xr, Double_t x[5]) const {
102 //
103 // This function returns external TOF track representation
104 //
105 xr=fX;
106 x[0]=GetY();
107 x[1]=GetZ();
108 x[2]=GetSnp();
109 x[3]=GetTgl();
110 x[4]=Get1Pt();
111}
112
113//_____________________________________________________________________________
114void AliTOFtrack::GetExternalCovariance(Double_t cc[15]) const {
115 //
116 // This function returns external representation of the covriance matrix.
117 //
118 Double_t a=GetConvConst();
119 Double_t c22=fX*fX*fCcc-2*fX*fCce+fCee;
120 Double_t c32=fX*fCct-fCte;
121 Double_t c20=fX*fCcy-fCey, c21=fX*fCcz-fCez, c42=fX*fCcc-fCce;
122
123 cc[0 ]=fCyy;
124 cc[1 ]=fCzy; cc[2 ]=fCzz;
125 cc[3 ]=c20; cc[4 ]=c21; cc[5 ]=c22;
126 cc[6 ]=fCty; cc[7 ]=fCtz; cc[8 ]=c32; cc[9 ]=fCtt;
127 cc[10]=fCcy*a; cc[11]=fCcz*a; cc[12]=c42*a; cc[13]=fCct*a; cc[14]=fCcc*a*a;
128
129}
130
131
132//_____________________________________________________________________________
133void AliTOFtrack::GetCovariance(Double_t cc[15]) const {
134
135 cc[0]=fCyy;
136 cc[1]=fCzy; cc[2]=fCzz;
137 cc[3]=fCey; cc[4]=fCez; cc[5]=fCee;
138 cc[6]=fCcy; cc[7]=fCcz; cc[8]=fCce; cc[9]=fCcc;
139 cc[10]=fCty; cc[11]=fCtz; cc[12]=fCte; cc[13]=fCct; cc[14]=fCtt;
140
141}
142
143
144//_____________________________________________________________________________
145Int_t AliTOFtrack::PropagateTo(Double_t xk,Double_t x0,Double_t rho)
146{
147 // Propagates a track of particle with mass=pm to a reference plane
148 // defined by x=xk through media of density=rho and radiationLength=x0
149
150 if (xk == fX) return 1;
151
152 if (TMath::Abs(fC*xk - fE) >= 0.90000) {
153 return 0;
154 }
155
156 // track Length measurement [SR, GSI, 17.02.2003]
157
158 Double_t oldX = fX, oldY = fY, oldZ = fZ;
159
160 Double_t x1=fX, x2=x1+(xk-x1), dx=x2-x1, y1=fY, z1=fZ;
161 Double_t c1=fC*x1 - fE;
162 if((c1*c1) > 1){
163 return 0;}
164 Double_t r1=sqrt(1.- c1*c1);
165 Double_t c2=fC*x2 - fE;
166 if((c2*c2) > 1) {
167 return 0;
168 }
169 Double_t r2=sqrt(1.- c2*c2);
170
171 fY += dx*(c1+c2)/(r1+r2);
172 fZ += dx*(c1+c2)/(c1*r2 + c2*r1)*fT;
173
174 //f = F - 1
175 Double_t rr=r1+r2, cc=c1+c2, xx=x1+x2;
176 Double_t f02=-dx*(2*rr + cc*(c1/r1 + c2/r2))/(rr*rr);
177 Double_t f04= dx*(rr*xx + cc*(c1*x1/r1+c2*x2/r2))/(rr*rr);
178 Double_t cr=c1*r2+c2*r1;
179 Double_t f12=-dx*fT*(2*cr + cc*(c2*c1/r1-r1 + c1*c2/r2-r2))/(cr*cr);
180 Double_t f13= dx*cc/cr;
181 Double_t f14=dx*fT*(cr*xx-cc*(r1*x2-c2*c1*x1/r1+r2*x1-c1*c2*x2/r2))/(cr*cr);
182
183 //b = C*ft
184 Double_t b00=f02*fCey + f04*fCcy, b01=f12*fCey + f14*fCcy + f13*fCty;
185 Double_t b10=f02*fCez + f04*fCcz, b11=f12*fCez + f14*fCcz + f13*fCtz;
186 Double_t b20=f02*fCee + f04*fCce, b21=f12*fCee + f14*fCce + f13*fCte;
187 Double_t b30=f02*fCte + f04*fCct, b31=f12*fCte + f14*fCct + f13*fCtt;
188 Double_t b40=f02*fCce + f04*fCcc, b41=f12*fCce + f14*fCcc + f13*fCct;
189
190 //a = f*b = f*C*ft
191 Double_t a00=f02*b20+f04*b40,a01=f02*b21+f04*b41,a11=f12*b21+f14*b41+f13*b31;
192
193 //F*C*Ft = C + (a + b + bt)
194 fCyy += a00 + 2*b00;
195 fCzy += a01 + b01 + b10;
196 fCey += b20;
197 fCty += b30;
198 fCcy += b40;
199 fCzz += a11 + 2*b11;
200 fCez += b21;
201 fCtz += b31;
202 fCcz += b41;
203
204 fX=x2;
205
206 //Multiple scattering ******************
207 Double_t d=sqrt((x1-fX)*(x1-fX)+(y1-fY)*(y1-fY)+(z1-fZ)*(z1-fZ));
208 Double_t p2=(1.+ GetTgl()*GetTgl())/(Get1Pt()*Get1Pt());
209 Double_t beta2=p2/(p2 + GetMass()*GetMass());
210 Double_t theta2=14.1*14.1/(beta2*p2*1e6)*d/x0*rho;
211
212 Double_t ey=fC*fX - fE, ez=fT;
213 Double_t xz=fC*ez, zz1=ez*ez+1, xy=fE+ey;
214
215 fCee += (2*ey*ez*ez*fE+1-ey*ey+ez*ez+fE*fE*ez*ez)*theta2;
216 fCte += ez*zz1*xy*theta2;
217 fCtt += zz1*zz1*theta2;
218 fCce += xz*ez*xy*theta2;
219 fCct += xz*zz1*theta2;
220 fCcc += xz*xz*theta2;
221 /*
222 Double_t dc22 = (1-ey*ey+xz*xz*fX*fX)*theta2;
223 Double_t dc32 = (xz*fX*zz1)*theta2;
224 Double_t dc33 = (zz1*zz1)*theta2;
225 Double_t dc42 = (xz*fX*xz)*theta2;
226 Double_t dc43 = (zz1*xz)*theta2;
227 Double_t dc44 = (xz*xz)*theta2;
228 fCee += dc22;
229 fCte += dc32;
230 fCtt += dc33;
231 fCce += dc42;
232 fCct += dc43;
233 fCcc += dc44;
234 */
235 //Energy losses************************
236 if((5940*beta2/(1-beta2+1e-10) - beta2) < 0){return 0;}
237
238 Double_t dE=0.153e-3/beta2*(log(5940*beta2/(1-beta2+1e-10)) - beta2)*d*rho;
239 if (x1 < x2) dE=-dE;
240 cc=fC;
241 fC*=(1.- sqrt(p2+GetMass()*GetMass())/p2*dE);
242 fE+=fX*(fC-cc);
243
244 // track time measurement [SR, GSI 17.02.2002]
245 if (x1 < x2)
246 if (IsStartedTimeIntegral()) {
247 Double_t l2 = (fX-oldX)*(fX-oldX) + (fY-oldY)*(fY-oldY) + (fZ-oldZ)*(fZ-oldZ);
248 AddTimeStep(TMath::Sqrt(l2));
249 }
250
251 return 1;
252}
253
254//_____________________________________________________________________________
255Int_t AliTOFtrack::PropagateToInnerTOF( Bool_t holes)
256{
257 // Propagates a track of particle with mass=pm to a reference plane
258 // defined by x=xk through media of density=rho and radiationLength=x0
259
260
261 Double_t ymax=AliTOFGeometry::RinTOF()*TMath::Tan(0.5*AliTOFGeometry::GetAlpha());
262 Bool_t skip = kFALSE;
263 Double_t y=this->GetYat(AliTOFGeometry::RinTOF(),skip);
264 if(skip){
265 return 0;
266 }
267 if (y > ymax) {
268 if (!this->Rotate(AliTOFGeometry::GetAlpha())) {
269 return 0;
270 }
271 } else if (y <-ymax) {
272 if (!this->Rotate(-AliTOFGeometry::GetAlpha())) {
273 return 0;
274 }
275 }
276
277
278 Double_t x = this->GetX();
279 Int_t nsteps=Int_t((370.-x)/0.5); // 0.5 cm Steps
280 for (Int_t istep=0;istep<nsteps;istep++){
281 Float_t xp = x+istep*0.5;
282 Double_t param[2];
283 this->GetPropagationParameters(holes,param);
284 this->PropagateTo(xp,param[0],param[1]);
285
286 }
287
288 if(!this->PropagateTo(AliTOFGeometry::RinTOF()))return 0;
289
290 return 1;
291
292}
293
294//_____________________________________________________________________________
295Int_t AliTOFtrack::Rotate(Double_t alpha)
296{
297 // Rotates track parameters in R*phi plane
298
299
300 fAlpha += alpha;
301 if (fAlpha<-TMath::Pi()) fAlpha += 2*TMath::Pi();
302 if (fAlpha>=TMath::Pi()) fAlpha -= 2*TMath::Pi();
303
304 Double_t x1=fX, y1=fY;
305 Double_t ca=cos(alpha), sa=sin(alpha);
306 Double_t r1=fC*fX - fE;
307
308 fX = x1*ca + y1*sa;
309 fY =-x1*sa + y1*ca;
310 if((r1*r1) > 1) return 0;
311 fE=fE*ca + (fC*y1 + sqrt(1.- r1*r1))*sa;
312
313 Double_t r2=fC*fX - fE;
314 if (TMath::Abs(r2) >= 0.90000) {
315 cerr<<" AliTOFtrack warning: Rotation failed !\n";
316 return 0;
317 }
318
319 if((r2*r2) > 1) return 0;
320 Double_t y0=fY + sqrt(1.- r2*r2)/fC;
321 if ((fY-y0)*fC >= 0.) {
322 cerr<<" AliTOFtrack warning: Rotation failed !!!\n";
323 return 0;
324 }
325
326 //f = F - 1
327 Double_t f00=ca-1, f24=(y1 - r1*x1/sqrt(1.- r1*r1))*sa,
328 f20=fC*sa, f22=(ca + sa*r1/sqrt(1.- r1*r1))-1;
329
330 //b = C*ft
331 Double_t b00=fCyy*f00, b02=fCyy*f20+fCcy*f24+fCey*f22;
332 Double_t b10=fCzy*f00, b12=fCzy*f20+fCcz*f24+fCez*f22;
333 Double_t b20=fCey*f00, b22=fCey*f20+fCce*f24+fCee*f22;
334 Double_t b30=fCty*f00, b32=fCty*f20+fCct*f24+fCte*f22;
335 Double_t b40=fCcy*f00, b42=fCcy*f20+fCcc*f24+fCce*f22;
336
337 //a = f*b = f*C*ft
338 Double_t a00=f00*b00, a02=f00*b02, a22=f20*b02+f24*b42+f22*b22;
339
340 //F*C*Ft = C + (a + b + bt)
341 fCyy += a00 + 2*b00;
342 fCzy += b10;
343 fCey += a02+b20+b02;
344 fCty += b30;
345 fCcy += b40;
346 fCez += b12;
347 fCte += b32;
348 fCee += a22 + 2*b22;
349 fCce += b42;
350
351 return 1;
352}
353
354
355
356//_________________________________________________________________________
357void AliTOFtrack::GetPxPyPz(Double_t& px, Double_t& py, Double_t& pz) const
358{
359 // Returns reconstructed track momentum in the global system.
360
361 Double_t pt=TMath::Abs(GetPt()); // GeV/c
362 Double_t r=fC*fX-fE;
363
364 Double_t y0;
365 if(r > 1) { py = pt; px = 0; }
366 else if(r < -1) { py = -pt; px = 0; }
367 else {
368 y0=fY + sqrt(1.- r*r)/fC;
369 px=-pt*(fY-y0)*fC; //cos(phi);
370 py=-pt*(fE-fX*fC); //sin(phi);
371 }
372 pz=pt*fT;
373 Double_t tmp=px*TMath::Cos(fAlpha) - py*TMath::Sin(fAlpha);
374 py=px*TMath::Sin(fAlpha) + py*TMath::Cos(fAlpha);
375 px=tmp;
376
377}
378
379//_________________________________________________________________________
380void AliTOFtrack::GetGlobalXYZ(Double_t& x, Double_t& y, Double_t& z) const
381{
382 // Returns reconstructed track coordinates in the global system.
383
384 x = fX; y = fY; z = fZ;
385 Double_t tmp=x*TMath::Cos(fAlpha) - y*TMath::Sin(fAlpha);
386 y=x*TMath::Sin(fAlpha) + y*TMath::Cos(fAlpha);
387 x=tmp;
388
389}
390
391//_________________________________________________________________________
392void AliTOFtrack::ResetCovariance() {
393 //
394 // Resets covariance matrix
395 //
396
397 fCyy*=10.;
398 fCzy=0.; fCzz*=10.;
399 fCey=0.; fCez=0.; fCee*=10.;
400 fCty=0.; fCtz=0.; fCte=0.; fCtt*=10.;
401 fCcy=0.; fCcz=0.; fCce=0.; fCct=0.; fCcc*=10.;
402}
403
404
405//_________________________________________________________________________
406void AliTOFtrack::ResetCovariance(Float_t mult) {
407 //
408 // Resets covariance matrix
409 //
410
411 fCyy*=mult;
412 fCzy*=0.; fCzz*=mult;
413 fCey*=0.; fCez*=0.; fCee*=mult;
414 fCty*=0.; fCtz*=0.; fCte*=0.; fCtt*=mult;
415 fCcy*=0.; fCcz*=0.; fCce*=0.; fCct*=0.; fCcc*=mult;
416}
417
418//_____________________________________________________________________________
419Int_t AliTOFtrack::Compare(const TObject *o) const {
420 //-----------------------------------------------------------------
421 // This function compares tracks according to the their curvature
422 //-----------------------------------------------------------------
423 AliTOFtrack *t=(AliTOFtrack*)o;
424 Double_t co=t->GetSigmaY2()*t->GetSigmaZ2();
425 Double_t c =GetSigmaY2()*GetSigmaZ2();
426 if (c>co) return 1;
427 else if (c<co) return -1;
428 return 0;
429}
430
431//_____________________________________________________________________________
432void AliTOFtrack::GetPropagationParameters(Bool_t holes, Double_t *param) {
433
434
435 Double_t thetamin = (90.-31.1) * TMath::Pi()/180.;
436 Double_t thetamax = (90.+31.1) * TMath::Pi()/180.;
437
438 Double_t zmin = -55.;
439 Double_t zmax = 55.;
440
441 Double_t rTPC = 261.53;
442 Double_t rTPCTRD = 294.5;
443 Double_t rTRD = 369.1;
444
445 Double_t x0TPC = 40.;
446 Double_t rhoTPC =0.06124;
447
448 Double_t x0Air = 36.66;
449 Double_t rhoAir =1.2931e-3;
450
451 Double_t x0TRD = 171.7;
452 Double_t rhoTRD =0.33;
453
454 Int_t isec = this->GetSector();
455 Double_t xtr,ytr,ztr;
456 this->GetGlobalXYZ(xtr,ytr,ztr);
457 Float_t thetatr = TMath::ATan2(TMath::Sqrt(xtr*xtr+ytr*ytr),ztr);
458
459 if(holes){
460 if (isec == 0 || isec == 1 || isec == 2 ) {
461 if( thetatr>=thetamin && thetatr<=thetamax){
462 x0TRD= x0Air;
463 rhoTRD = rhoAir;
464 }
465 }
466 if (isec == 11 || isec == 12 || isec == 13 || isec == 14 || isec == 15 ) {
467 if( ztr>=zmin && ztr<=zmax){
468 x0TRD= x0Air;
469 rhoTRD = rhoAir;
470 }
471 }
472 }
473
474 if(this->GetX() <= rTPC)
475 {param[0]=x0TPC;param[1]=rhoTPC;}
476 else if(this->GetX() > rTPC && this->GetX() < rTPCTRD)
477 {param[0]=x0Air;param[1]=rhoAir;}
478 else if(this->GetX() >= rTPCTRD && this->GetX() < rTRD)
479 {param[0]=x0TRD;param[1]=rhoTRD;}
480 else if(this->GetX() >= rTRD )
481 {param[0]=x0Air;param[1]=rhoAir;}
482
483}