]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TRD/AliTRDtrack.cxx
J.Chudoba's changes merged correctly with the HEAD
[u/mrichter/AliRoot.git] / TRD / AliTRDtrack.cxx
CommitLineData
46d29e70 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$
3ab6f951 18Revision 1.7 2001/05/28 17:07:58 hristov
19Last minute changes; ExB correction in AliTRDclusterizerV1; taking into account of material in G10 TEC frames and material between TEC planes (C.Blume,S.Sedykh)
20
2a941f4e 21Revision 1.4 2000/12/08 16:07:02 cblume
22Update of the tracking by Sergei
23
bbf92647 24Revision 1.3 2000/10/15 23:40:01 cblume
25Remove AliTRDconst
26
0e9c2ad5 27Revision 1.2 2000/10/06 16:49:46 cblume
28Made Getters const
29
46d29e70 30Revision 1.1.2.1 2000/09/22 14:47:52 cblume
31Add the tracking code
32
33*/
34
35#include <iostream.h>
36
37#include <TObject.h>
38
39#include "AliTRD.h"
46d29e70 40#include "AliTRDgeometry.h"
41#include "AliTRDcluster.h"
42#include "AliTRDtrack.h"
43
44ClassImp(AliTRDtrack)
45
46
47//_____________________________________________________________________________
48
a819a5f7 49AliTRDtrack::AliTRDtrack(const AliTRDcluster *c, UInt_t index,
50const Double_t xx[5], const Double_t cc[15], Double_t xref, Double_t alpha) {
46d29e70 51 //-----------------------------------------------------------------
52 // This is the main track constructor.
53 //-----------------------------------------------------------------
54 fLab=-1;
55 fChi2=0.;
56 fdEdx=0.;
57
58 fAlpha=alpha;
59 fX=xref;
60
61 fY=xx[0]; fZ=xx[1]; fC=xx[2]; fE=xx[3]; fT=xx[4];
62
63 fCyy=cc[0];
64 fCzy=cc[1]; fCzz=cc[2];
65 fCcy=cc[3]; fCcz=cc[4]; fCcc=cc[5];
66 fCey=cc[6]; fCez=cc[7]; fCec=cc[8]; fCee=cc[9];
67 fCty=cc[10]; fCtz=cc[11]; fCtc=cc[12]; fCte=cc[13]; fCtt=cc[14];
68
69 fN=0;
a819a5f7 70 fIndex[fN]=index;
a2b90f83 71
a819a5f7 72 Float_t q = c->GetQ();
73 Double_t s = fX*fC - fE, t=fT;
74 q *= TMath::Sqrt((1-s*s)/(1+t*t));
75 fdQdl[fN++] = q;
46d29e70 76}
77
78//_____________________________________________________________________________
79AliTRDtrack::AliTRDtrack(const AliTRDtrack& t) {
80 //
81 // Copy constructor.
82 //
83
84 fLab=t.fLab;
85
86 fChi2=t.fChi2;
87 fdEdx=t.fdEdx;
88
89 fAlpha=t.fAlpha;
90 fX=t.fX;
91
92 fY=t.fY; fZ=t.fZ; fC=t.fC; fE=t.fE; fT=t.fT;
93
94 fCyy=t.fCyy;
95 fCzy=t.fCzy; fCzz=t.fCzz;
96 fCcy=t.fCcy; fCcz=t.fCcz; fCcc=t.fCcc;
97 fCey=t.fCey; fCez=t.fCez; fCec=t.fCec; fCee=t.fCee;
98 fCty=t.fCty; fCtz=t.fCtz; fCtc=t.fCtc; fCte=t.fCte; fCtt=t.fCtt;
99
100 fN=t.fN;
a819a5f7 101 for (Int_t i=0; i<fN; i++) {
102 fIndex[i]=t.fIndex[i];
103 fdQdl[i]=t.fdQdl[i];
104 }
46d29e70 105}
106
107//_____________________________________________________________________________
108void AliTRDtrack::GetCovariance(Double_t cc[15]) const {
109 cc[0]=fCyy;
110 cc[1]=fCzy; cc[2]=fCzz;
111 cc[3]=fCcy; cc[4]=fCcz; cc[5]=fCcc;
112 cc[6]=fCey; cc[7]=fCez; cc[8]=fCec; cc[9]=fCee;
113 cc[10]=fCty; cc[11]=fCtz; cc[12]=fCtc; cc[13]=fCte; cc[14]=fCtt;
114}
115
116//_____________________________________________________________________________
2a941f4e 117Int_t AliTRDtrack::Compare(const TObject *o) const {
46d29e70 118
119// Compares tracks according to their Y2
120
121 AliTRDtrack *t=(AliTRDtrack*)o;
122 // Double_t co=t->GetSigmaY2();
123 // Double_t c =GetSigmaY2();
124
125 Double_t co=TMath::Abs(t->GetC());
126 Double_t c =TMath::Abs(GetC());
127
128 if (c>co) return 1;
129 else if (c<co) return -1;
130 return 0;
131}
132
a819a5f7 133//_____________________________________________________________________________
134void AliTRDtrack::CookdEdx(Double_t low, Double_t up) {
135 //-----------------------------------------------------------------
136 // Calculates dE/dX within the "low" and "up" cuts.
137 //-----------------------------------------------------------------
138 Int_t i;
139 Int_t nc=GetNclusters();
140
141 Float_t sorted[200];
142 for (i=0; i<200; i++) sorted[i]=fdQdl[i];
143
144 Int_t swap;
145 do {
146 swap=0;
147 for (i=0; i<nc-1; i++) {
148 if (sorted[i]<=sorted[i+1]) continue;
149 Float_t tmp=sorted[i];
150 sorted[i]=sorted[i+1]; sorted[i+1]=tmp;
151 swap++;
152 }
153 } while (swap);
154
155 Int_t nl=Int_t(low*nc), nu=Int_t(up*nc);
156 Float_t dedx=0;
157 for (i=nl; i<=nu; i++) dedx += sorted[i];
158 dedx /= (nu-nl+1);
159 SetdEdx(dedx);
160}
161
162
163
46d29e70 164//_____________________________________________________________________________
165Int_t AliTRDtrack::PropagateTo(Double_t xk,Double_t x0,Double_t rho,Double_t pm)
166{
167 // Propagates a track of particle with mass=pm to a reference plane
168 // defined by x=xk through media of density=rho and radiationLength=x0
169
170 if (TMath::Abs(fC*xk - fE) >= 0.99999) {
171 if (fN>4) cerr<<fN<<" AliTRDtrack warning: Propagation failed !\n";
172 return 0;
173 }
174
175 Double_t x1=fX, x2=x1+(xk-x1), dx=x2-x1, y1=fY, z1=fZ;
176 Double_t c1=fC*x1 - fE, r1=sqrt(1.- c1*c1);
177 Double_t c2=fC*x2 - fE, r2=sqrt(1.- c2*c2);
178
179 fY += dx*(c1+c2)/(r1+r2);
180 fZ += dx*(c1+c2)/(c1*r2 + c2*r1)*fT;
181
182 //f = F - 1
183 Double_t rr=r1+r2, cc=c1+c2, xx=x1+x2;
184 Double_t f02= dx*(rr*xx + cc*(c1*x1/r1+c2*x2/r2))/(rr*rr);
185 Double_t f03=-dx*(2*rr + cc*(c1/r1 + c2/r2))/(rr*rr);
186 Double_t cr=c1*r2+c2*r1;
187 Double_t f12= dx*fT*(cr*xx-cc*(r1*x2-c2*c1*x1/r1+r2*x1-c1*c2*x2/r2))/(cr*cr);
188 Double_t f13=-dx*fT*(2*cr + cc*(c2*c1/r1-r1 + c1*c2/r2-r2))/(cr*cr);
189 Double_t f14= dx*cc/cr;
190
191 //b = C*ft
192 Double_t b00=f02*fCcy + f03*fCey, b01=f12*fCcy + f13*fCey + f14*fCty;
193 Double_t b10=f02*fCcz + f03*fCez, b11=f12*fCcz + f13*fCez + f14*fCtz;
194 Double_t b20=f02*fCcc + f03*fCec, b21=f12*fCcc + f13*fCec + f14*fCtc;
195 Double_t b30=f02*fCec + f03*fCee, b31=f12*fCec + f13*fCee + f14*fCte;
196 Double_t b40=f02*fCtc + f03*fCte, b41=f12*fCtc + f13*fCte + f14*fCtt;
197
198 //a = f*b = f*C*ft
199 Double_t a00=f02*b20+f03*b30,a01=f02*b21+f03*b31,a11=f12*b21+f13*b31+f14*b41;
200
201 //F*C*Ft = C + (a + b + bt)
202 fCyy += a00 + 2*b00;
203 fCzy += a01 + b01 + b10;
204 fCcy += b20;
205 fCey += b30;
206 fCty += b40;
207 fCzz += a11 + 2*b11;
208 fCcz += b21;
209 fCez += b31;
210 fCtz += b41;
211
212 fX=x2;
213
214
215 //Multiple scattering ******************
216
217 Double_t d=sqrt((x1-fX)*(x1-fX)+(y1-fY)*(y1-fY)+(z1-fZ)*(z1-fZ));
218 Double_t p2=GetPt()*GetPt()*(1.+fT*fT);
219 Double_t beta2=p2/(p2 + pm*pm);
220
221 Double_t ey=fC*fX - fE, ez=fT;
222 Double_t xz=fC*ez, zz1=ez*ez+1, xy=fE+ey;
223
224 Double_t theta2=14.1*14.1/(beta2*p2*1e6)*d/x0*rho;
225 fCcc += xz*xz*theta2;
226 fCec += xz*ez*xy*theta2;
227 fCtc += xz*zz1*theta2;
228 fCee += (2*ey*ez*ez*fE+1-ey*ey+ez*ez+fE*fE*ez*ez)*theta2;
229 fCte += ez*zz1*xy*theta2;
230 fCtt += zz1*zz1*theta2;
231
232
233 //Energy losses************************
234
235 Double_t dE=0.153e-3/beta2*(log(5940*beta2/(1-beta2)) - beta2)*d*rho;
236 if (x1 < x2) dE=-dE;
237 fC*=(1.- sqrt(p2+pm*pm)/p2*dE);
238 //fE*=(1.- sqrt(p2+pm*pm)/p2*dE);
239
240 return 1;
241
242}
243
244
46d29e70 245//_____________________________________________________________________________
246void AliTRDtrack::Update(const AliTRDcluster *c, Double_t chisq, UInt_t index)
247{
248 // Assignes found cluster to the track and updates track information
249
a819a5f7 250 Double_t r00=c->GetSigmaY2(), r01=0., r11=c->GetSigmaZ2();
46d29e70 251 r00+=fCyy; r01+=fCzy; r11+=fCzz;
252 Double_t det=r00*r11 - r01*r01;
253 Double_t tmp=r00; r00=r11/det; r11=tmp/det; r01=-r01/det;
254
255 Double_t k00=fCyy*r00+fCzy*r01, k01=fCyy*r01+fCzy*r11;
256 Double_t k10=fCzy*r00+fCzz*r01, k11=fCzy*r01+fCzz*r11;
257 Double_t k20=fCcy*r00+fCcz*r01, k21=fCcy*r01+fCcz*r11;
258 Double_t k30=fCey*r00+fCez*r01, k31=fCey*r01+fCez*r11;
259 Double_t k40=fCty*r00+fCtz*r01, k41=fCty*r01+fCtz*r11;
260
261 Double_t dy=c->GetY() - fY, dz=c->GetZ() - fZ;
262 Double_t cur=fC + k20*dy + k21*dz, eta=fE + k30*dy + k31*dz;
263 if (TMath::Abs(cur*fX-eta) >= 0.99999) {
264 if (fN>4) cerr<<fN<<" AliTRDtrack warning: Filtering failed !\n";
265 return;
266 }
267
268 fY += k00*dy + k01*dz;
269 fZ += k10*dy + k11*dz;
270 fC = cur;
271 fE = eta;
272 fT += k40*dy + k41*dz;
273
274 Double_t c01=fCzy, c02=fCcy, c03=fCey, c04=fCty;
275 Double_t c12=fCcz, c13=fCez, c14=fCtz;
276
277 fCyy-=k00*fCyy+k01*fCzy; fCzy-=k00*c01+k01*fCzz;
278 fCcy-=k00*c02+k01*c12; fCey-=k00*c03+k01*c13;
279 fCty-=k00*c04+k01*c14;
280
281 fCzz-=k10*c01+k11*fCzz;
282 fCcz-=k10*c02+k11*c12; fCez-=k10*c03+k11*c13;
283 fCtz-=k10*c04+k11*c14;
284
285 fCcc-=k20*c02+k21*c12; fCec-=k20*c03+k21*c13;
286 fCtc-=k20*c04+k21*c14;
287
288 fCee-=k30*c03+k31*c13;
289 fCte-=k30*c04+k31*c14;
290
291 fCtt-=k40*c04+k41*c14;
292
a819a5f7 293 fIndex[fN]=index;
294
295 Float_t q = c->GetQ();
296 Double_t s = fX*fC - fE, t=fT;
297 q *= TMath::Sqrt((1-s*s)/(1+t*t));
298 fdQdl[fN++] = q;
299
46d29e70 300 fChi2 += chisq;
301
302 // cerr<<"in update: fIndex["<<fN<<"] = "<<index<<endl;
303}
304
305//_____________________________________________________________________________
306Int_t AliTRDtrack::Rotate(Double_t alpha)
307{
308 // Rotates track parameters in R*phi plane
309
310 fAlpha += alpha;
311
312 Double_t x1=fX, y1=fY;
313 Double_t ca=cos(alpha), sa=sin(alpha);
314 Double_t r1=fC*fX - fE;
315
316 fX = x1*ca + y1*sa;
317 fY=-x1*sa + y1*ca;
318 fE=fE*ca + (fC*y1 + sqrt(1.- r1*r1))*sa;
319
320 Double_t r2=fC*fX - fE;
321 if (TMath::Abs(r2) >= 0.99999) {
322 if (fN>4) cerr<<fN<<" AliTRDtrack warning: Rotation failed !\n";
323 return 0;
324 }
325
326 Double_t y0=fY + sqrt(1.- r2*r2)/fC;
327 if ((fY-y0)*fC >= 0.) {
328 if (fN>4) cerr<<fN<<" AliTRDtrack warning: Rotation failed !!!\n";
329 return 0;
330 }
331
332 //f = F - 1
333 Double_t f00=ca-1, f32=(y1 - r1*x1/sqrt(1.- r1*r1))*sa,
334 f30=fC*sa, f33=(ca + sa*r1/sqrt(1.- r1*r1))-1;
335
336 //b = C*ft
337 Double_t b00=fCyy*f00, b03=fCyy*f30+fCcy*f32+fCey*f33;
338 Double_t b10=fCzy*f00, b13=fCzy*f30+fCcz*f32+fCez*f33;
339 Double_t b20=fCcy*f00, b23=fCcy*f30+fCcc*f32+fCec*f33;
340 Double_t b30=fCey*f00, b33=fCey*f30+fCec*f32+fCee*f33;
341 Double_t b40=fCty*f00, b43=fCty*f30+fCtc*f32+fCte*f33;
342
343 //a = f*b = f*C*ft
344 Double_t a00=f00*b00, a03=f00*b03, a33=f30*b03+f32*b23+f33*b33;
345
346 // *** Double_t dy2=fCyy;
347
348 //F*C*Ft = C + (a + b + bt)
349 fCyy += a00 + 2*b00;
350 fCzy += b10;
351 fCcy += b20;
352 fCey += a03+b30+b03;
353 fCty += b40;
354 fCez += b13;
355 fCec += b23;
356 fCee += a33 + 2*b33;
357 fCte += b43;
358
359 // *** fCyy+=dy2*sa*sa*r1*r1/(1.- r1*r1);
360 // *** fCzz+=d2y*sa*sa*fT*fT/(1.- r1*r1);
361
362 return 1;
363}
364
365
46d29e70 366//_____________________________________________________________________________
367Double_t AliTRDtrack::GetPredictedChi2(const AliTRDcluster *c) const
368{
a819a5f7 369 /*
370 Double_t r00=c->GetSigmaY2(), r01=0., r11=c->GetSigmaZ2();
46d29e70 371 r00+=fCyy; r01+=fCzy; r11+=fCzz;
372
373 Double_t det=r00*r11 - r01*r01;
374 if (TMath::Abs(det) < 1.e-10) {
375 if (fN>4) cerr<<fN<<" AliTRDtrack warning: Singular matrix !\n";
376 return 1e10;
377 }
378 Double_t tmp=r00; r00=r11; r11=tmp; r01=-r01;
379
380 Double_t dy=c->GetY() - fY, dz=c->GetZ() - fZ;
381
382 return (dy*r00*dy + 2*r01*dy*dz + dz*r11*dz)/det;
a819a5f7 383 */
384
385 Double_t dy=c->GetY() - fY;
386 Double_t r00=c->GetSigmaY2();
387
388 return (dy*dy)/r00;
389
46d29e70 390}
391
392
393//_________________________________________________________________________
394void AliTRDtrack::GetPxPyPz(Double_t& px, Double_t& py, Double_t& pz) const
395{
396 // Returns reconstructed track momentum in the global system.
397
398 Double_t pt=TMath::Abs(GetPt()); // GeV/c
399 Double_t r=fC*fX-fE;
400 Double_t y0=fY + sqrt(1.- r*r)/fC;
401 px=-pt*(fY-y0)*fC; //cos(phi);
402 py=-pt*(fE-fX*fC); //sin(phi);
403 pz=pt*fT;
404 Double_t tmp=px*TMath::Cos(fAlpha) - py*TMath::Sin(fAlpha);
405 py=px*TMath::Sin(fAlpha) + py*TMath::Cos(fAlpha);
406 px=tmp;
407
408}
409
410//____________________________________________________________________________
411void AliTRDtrack::Streamer(TBuffer &R__b)
412{
3ab6f951 413 Int_t i;
414
46d29e70 415 if (R__b.IsReading()) {
416 Version_t R__v = R__b.ReadVersion(); if (R__v) { }
417 TObject::Streamer(R__b);
418 R__b >> fLab;
419 R__b >> fChi2;
420 R__b >> fdEdx;
421 R__b >> fAlpha;
422 R__b >> fX;
423 R__b >> fY;
424 R__b >> fZ;
425 R__b >> fC;
426 R__b >> fE;
427 R__b >> fT;
428 R__b >> fCyy;
429 R__b >> fCzy;
430 R__b >> fCzz;
431 R__b >> fCcy;
432 R__b >> fCcz;
433 R__b >> fCcc;
434 R__b >> fCey;
435 R__b >> fCez;
436 R__b >> fCec;
437 R__b >> fCee;
438 R__b >> fCty;
439 R__b >> fCtz;
440 R__b >> fCtc;
441 R__b >> fCte;
442 R__b >> fCtt;
443 R__b >> fN;
3ab6f951 444 for (i=0; i<fN; i++) R__b >> fIndex[i];
445 for (i=0; i<fN; i++) R__b >> fdQdl[i];
46d29e70 446 } else {
447 R__b.WriteVersion(AliTRDtrack::IsA());
448 TObject::Streamer(R__b);
449 R__b << fLab;
450 R__b << fChi2;
451 R__b << fdEdx;
452 R__b << fAlpha;
453 R__b << fX;
454 R__b << fY;
455 R__b << fZ;
456 R__b << fC;
457 R__b << fE;
458 R__b << fT;
459 R__b << fCyy;
460 R__b << fCzy;
461 R__b << fCzz;
462 R__b << fCcy;
463 R__b << fCcz;
464 R__b << fCcc;
465 R__b << fCey;
466 R__b << fCez;
467 R__b << fCec;
468 R__b << fCee;
469 R__b << fCty;
470 R__b << fCtz;
471 R__b << fCtc;
472 R__b << fCte;
473 R__b << fCtt;
474 R__b << fN;
3ab6f951 475 for (i=0; i<fN; i++) R__b << fIndex[i];
476 for (i=0; i<fN; i++) R__b << fdQdl[i];
46d29e70 477 }
478}
479
46d29e70 480