]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliESDtrack.cxx
Corrections to obey the coding conventions
[u/mrichter/AliRoot.git] / STEER / AliESDtrack.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 ESD track class
18 //   ESD = Event Summary Data
19 //   This is the class to deal with during the phisical analysis of data
20 //      Origin: Iouri Belikov, CERN
21 //      e-mail: Jouri.Belikov@cern.ch
22 //-----------------------------------------------------------------
23
24 #include "TMath.h"
25
26 #include "AliESDtrack.h"
27 #include "AliKalmanTrack.h"
28 #include "../ITS/AliITStrackV2.h"
29
30 ClassImp(AliESDtrack)
31
32 //_______________________________________________________________________
33 AliESDtrack::AliESDtrack() : 
34 fFlags(0),
35 fLabel(0),
36 fTrackLength(0),
37 fStopVertex(0),
38 fRalpha(0),
39 fRx(0),
40 fITSchi2(0),
41 fITSncls(0),
42 fITSsignal(0),
43 fVertexX(0),
44 fVertexY(0),
45 fVertexZ(0),
46 fVertexPx(0),
47 fVertexPy(0),
48 fVertexPz(0),
49 fVertex(kFALSE),
50 fTPCchi2(0),
51 fTPCncls(0),
52 fTPCsignal(0),
53 fTRDchi2(0),
54 fTRDncls(0),
55 fTRDsignal(0),
56 fTOFchi2(0),
57 fTOFindex(0),
58 fTOFsignal(0)
59 {
60   //
61   // The default ESD constructor 
62   //
63   for (Int_t i=0; i<kSPECIES; i++) {
64     fTrackTime[i]=0;
65     fR[i]=0;
66     fITSr[i]=0;
67     fTPCr[i]=0;
68     fTRDr[i]=0;
69     fTOFr[i]=0;
70   }
71   for (Int_t i=0; i<5; fRp[i++]);
72   for (Int_t i=0; i<15; fRc[i++]);
73   for (Int_t i=0; i<6; fITSindex[i++]);
74   for (Int_t i=0; i<180; fTPCindex[i++]);
75 }
76
77 //_______________________________________________________________________
78 Float_t AliESDtrack::GetMass() const {
79   // Returns the mass of the most probable particle type
80   Float_t max=0.;
81   Int_t k=-1;
82   for (Int_t i=0; i<kSPECIES; i++) {
83     if (fR[i]>max) {k=i; max=fR[i];}
84   }
85   if (k==0) return 0.00051;
86   if (k==1) return 0.10566;
87   if (k==2||k==-1) return 0.13957;
88   if (k==3) return 0.49368;
89   if (k==4) return 0.93827;
90   Warning("GetMass()","Undefined mass !");
91   return 0.13957;
92 }
93
94 //_______________________________________________________________________
95 Bool_t AliESDtrack::UpdateTrackParams(AliKalmanTrack *t, ULong_t flags) {
96   //
97   // This function updates track's running parameters 
98   //
99   switch (flags) {
100     
101   case kITSin:
102   case kITSout: 
103   case kITSrefit:
104     fITSncls=t->GetNumberOfClusters();
105     fITSchi2=t->GetChi2();
106     for (Int_t i=0;i<fITSncls;i++) fITSindex[i]=t->GetClusterIndex(i);
107     fITSsignal=t->GetPIDsignal();
108     break;
109     
110   case kTPCin: case kTPCout: case kTPCrefit:
111     fTPCncls=t->GetNumberOfClusters();
112     fTPCchi2=t->GetChi2();
113     for (Int_t i=0;i<fTPCncls;i++) fTPCindex[i]=t->GetClusterIndex(i);
114     fTPCsignal=t->GetPIDsignal();
115     {Double_t mass=t->GetMass();    // preliminary mass setting 
116     if (mass>0.5) fR[4]=1.;         //        used by
117     else if (mass<0.4) fR[2]=1.;    // the ITS reconstruction
118     else fR[3]=1.;}                 //
119     break;
120   case kTRDin: case kTRDout: case kTRDrefit:
121     fTRDncls=t->GetNumberOfClusters();
122     fTRDchi2=t->GetChi2();
123     fTRDsignal=t->GetPIDsignal();
124     break;
125   default: 
126     Error("UpdateTrackParams()","Wrong flag !\n");
127     return kFALSE;
128   }
129
130   SetStatus(flags);
131   fLabel=t->GetLabel();
132
133   if (t->IsStartedTimeIntegral()) {
134     SetStatus(kTIME);
135     Double_t times[10];t->GetIntegratedTimes(times); SetIntegratedTimes(times);
136     SetIntegratedLength(t->GetIntegratedLength());
137   }
138
139   fRalpha=t->GetAlpha();
140   t->GetExternalParameters(fRx,fRp);
141   t->GetExternalCovariance(fRc);
142   
143   if (flags == kITSin)
144    {
145      AliITStrackV2* itstrack = dynamic_cast<AliITStrackV2*>(t);
146      if (itstrack)
147       {
148         itstrack->PropagateTo(3.,0.0028,65.19);
149         itstrack->PropagateToVertex();
150         
151         Double_t ralpha=t->GetAlpha();
152         Double_t rx;      // X-coordinate of the track reference plane 
153         Double_t rp[5];   // external track parameters  
154         t->GetExternalParameters(rx,rp);
155    
156         Double_t phi=TMath::ASin(rp[2]) + ralpha;
157         Double_t pt=1./TMath::Abs(rp[4]);
158         Double_t r=TMath::Sqrt(rx*rx + rp[0]*rp[0]);
159         
160         fVertexX=r*TMath::Cos(phi); 
161         fVertexY=r*TMath::Sin(phi); 
162         fVertexZ=rp[1]; 
163         
164         fVertexPx = pt*TMath::Cos(phi); 
165         fVertexPy = pt*TMath::Sin(phi); 
166         fVertexPz = pt*rp[3]; 
167         fVertex = kTRUE;
168       }
169    }
170   
171   return kTRUE;
172 }
173
174 //_______________________________________________________________________
175 void AliESDtrack::GetExternalParametersAt(Double_t x, Double_t p[5]) const {
176   //---------------------------------------------------------------------
177   // This function returns external representation of the track parameters
178   // at the plane x
179   //---------------------------------------------------------------------
180   Double_t dx=x-fRx;
181   Double_t c=fRp[4]/AliKalmanTrack::GetConvConst();
182   Double_t f1=fRp[2], f2=f1 + c*dx;
183   Double_t r1=sqrt(1.- f1*f1), r2=sqrt(1.- f2*f2);    
184
185   p[0]=fRp[0]+dx*(f1+f2)/(r1+r2);
186   p[1]=fRp[1]+dx*(f1+f2)/(f1*r2 + f2*r1)*fRp[3];
187   p[2]=fRp[2]+dx*c;
188   p[3]=fRp[3];
189   p[4]=fRp[4];
190 }
191
192 //_______________________________________________________________________
193 void AliESDtrack::GetExternalParameters(Double_t &x, Double_t p[5]) const {
194   //---------------------------------------------------------------------
195   // This function returns external representation of the track parameters
196   //---------------------------------------------------------------------
197   x=fRx;
198   for (Int_t i=0; i<5; i++) p[i]=fRp[i];
199 }
200
201 Double_t AliESDtrack::GetP() const {
202   //---------------------------------------------------------------------
203   // This function returns the track momentum
204   //---------------------------------------------------------------------
205   Double_t lam=TMath::ATan(fRp[3]);
206   Double_t pt=1./TMath::Abs(fRp[4]);
207   return pt/TMath::Cos(lam);
208 }
209
210 void AliESDtrack::GetPxPyPz(Double_t *p) const {
211   //---------------------------------------------------------------------
212   // This function returns the global track momentum components
213   //---------------------------------------------------------------------
214   Double_t phi=TMath::ASin(fRp[2]) + fRalpha;
215   Double_t pt=1./TMath::Abs(fRp[4]);
216   p[0]=pt*TMath::Cos(phi); p[1]=pt*TMath::Sin(phi); p[2]=pt*fRp[3]; 
217 }
218
219 void AliESDtrack::GetXYZ(Double_t *xyz) const {
220   //---------------------------------------------------------------------
221   // This function returns the global track position
222   //---------------------------------------------------------------------
223   Double_t phi=TMath::ASin(fRp[2]) + fRalpha;
224   Double_t r=TMath::Sqrt(fRx*fRx + fRp[0]*fRp[0]);
225   xyz[0]=r*TMath::Cos(phi); xyz[1]=r*TMath::Sin(phi); xyz[2]=fRp[1]; 
226 }
227
228 //_______________________________________________________________________
229 void AliESDtrack::GetExternalCovariance(Double_t c[15]) const {
230   //---------------------------------------------------------------------
231   // This function returns external representation of the cov. matrix
232   //---------------------------------------------------------------------
233   for (Int_t i=0; i<15; i++) c[i]=fRc[i];
234 }
235
236 //_______________________________________________________________________
237 void AliESDtrack::GetIntegratedTimes(Double_t *times) const {
238   // Returns the array with integrated times for each particle hypothesis
239   for (Int_t i=0; i<kSPECIES; i++) times[i]=fTrackTime[i];
240 }
241
242 //_______________________________________________________________________
243 void AliESDtrack::SetIntegratedTimes(const Double_t *times) {
244   // Sets the array with integrated times for each particle hypotesis
245   for (Int_t i=0; i<kSPECIES; i++) fTrackTime[i]=times[i];
246 }
247
248 //_______________________________________________________________________
249 void AliESDtrack::SetITSpid(const Double_t *p) {
250   // Sets values for the probability of each particle type (in ITS)
251   for (Int_t i=0; i<kSPECIES; i++) fITSr[i]=p[i];
252   SetStatus(AliESDtrack::kITSpid);
253 }
254
255 //_______________________________________________________________________
256 void AliESDtrack::GetITSpid(Double_t *p) const {
257   // Gets the probability of each particle type (in ITS)
258   for (Int_t i=0; i<kSPECIES; i++) p[i]=fITSr[i];
259 }
260
261 //_______________________________________________________________________
262 Int_t AliESDtrack::GetITSclusters(UInt_t *idx) const {
263   //---------------------------------------------------------------------
264   // This function returns indices of the assgined ITS clusters 
265   //---------------------------------------------------------------------
266   for (Int_t i=0; i<fITSncls; i++) idx[i]=fITSindex[i];
267   return fITSncls;
268 }
269
270 //_______________________________________________________________________
271 Int_t AliESDtrack::GetTPCclusters(Int_t *idx) const {
272   //---------------------------------------------------------------------
273   // This function returns indices of the assgined ITS clusters 
274   //---------------------------------------------------------------------
275   for (Int_t i=0; i<180; i++) idx[i]=fTPCindex[i];  // MI I prefer some constant
276   return fTPCncls;
277 }
278
279 //_______________________________________________________________________
280 void AliESDtrack::SetTPCpid(const Double_t *p) {  
281   // Sets values for the probability of each particle type (in TPC)
282   for (Int_t i=0; i<kSPECIES; i++) fTPCr[i]=p[i];
283   SetStatus(AliESDtrack::kTPCpid);
284 }
285
286 //_______________________________________________________________________
287 void AliESDtrack::GetTPCpid(Double_t *p) const {
288   // Gets the probability of each particle type (in TPC)
289   for (Int_t i=0; i<kSPECIES; i++) p[i]=fTPCr[i];
290 }
291
292 //_______________________________________________________________________
293 void AliESDtrack::SetTRDpid(const Double_t *p) {  
294   // Sets values for the probability of each particle type (in TRD)
295   for (Int_t i=0; i<kSPECIES; i++) fTRDr[i]=p[i];
296   SetStatus(AliESDtrack::kTRDpid);
297 }
298
299 //_______________________________________________________________________
300 void AliESDtrack::GetTRDpid(Double_t *p) const {
301   // Gets the probability of each particle type (in TRD)
302   for (Int_t i=0; i<kSPECIES; i++) p[i]=fTRDr[i];
303 }
304
305 //_______________________________________________________________________
306 void    AliESDtrack::SetTRDpid(Int_t iSpecies, Float_t p)
307 {
308   // Sets the probability of particle type iSpecies to p (in TRD)
309   fTRDr[iSpecies] = p;
310 }
311
312 Float_t AliESDtrack::GetTRDpid(Int_t iSpecies) const
313 {
314   // Returns the probability of particle type iSpecies (in TRD)
315   return fTRDr[iSpecies];
316 }
317
318 //_______________________________________________________________________
319 void AliESDtrack::SetTOFpid(const Double_t *p) {  
320   // Sets the probability of each particle type (in TOF)
321   for (Int_t i=0; i<kSPECIES; i++) fTOFr[i]=p[i];
322   SetStatus(AliESDtrack::kTOFpid);
323 }
324
325 //_______________________________________________________________________
326 void AliESDtrack::GetTOFpid(Double_t *p) const {
327   // Gets probabilities of each particle type (in TOF)
328   for (Int_t i=0; i<kSPECIES; i++) p[i]=fTOFr[i];
329 }
330
331 //_______________________________________________________________________
332 void AliESDtrack::SetESDpid(const Double_t *p) {  
333   // Sets the probability of each particle type for the ESD track
334   for (Int_t i=0; i<kSPECIES; i++) fR[i]=p[i];
335   SetStatus(AliESDtrack::kESDpid);
336 }
337
338 //_______________________________________________________________________
339 void AliESDtrack::GetESDpid(Double_t *p) const {
340   // Gets probability of each particle type for the ESD track
341   for (Int_t i=0; i<kSPECIES; i++) p[i]=fR[i];
342 }
343
344 void AliESDtrack::GetVertexXYZ(Double_t& x,Double_t& y, Double_t&z) const
345 {
346 //returns track position in DCA to vertex  
347   x = fVertexX;
348   y = fVertexY;
349   z = fVertexZ;
350 }
351 void AliESDtrack::GetVertexPxPyPz(Double_t& px,Double_t& py, Double_t& pz) const
352 {
353 //returns track momentum in DCA to vertex  
354   px = fVertexPx;
355   py = fVertexPy;
356   pz = fVertexPz;
357 }