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