]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliESDtrack.cxx
TRD included in the ESD chain (Yu.Belikov)
[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
29 ClassImp(AliESDtrack)
30
31 //_______________________________________________________________________
32 AliESDtrack::AliESDtrack() : 
33 fFlags(0),
34 fLabel(0),
35 fTrackLength(0),
36 fStopVertex(0),
37 fRalpha(0),
38 fRx(0),
39 fITSchi2(0),
40 fITSncls(0),
41 fITSsignal(0),
42 fVertexX(0),
43 fVertexY(0),
44 fVertexZ(0),
45 fVertexPx(0),
46 fVertexPy(0),
47 fVertexPz(0),
48 fVertex(kFALSE),
49 fTPCchi2(0),
50 fTPCncls(0),
51 fTPCsignal(0),
52 fTRDchi2(0),
53 fTRDncls(0),
54 fTRDsignal(0),
55 fTOFchi2(0),
56 fTOFindex(0),
57 fTOFsignal(-1)
58 {
59   //
60   // The default ESD constructor 
61   //
62   for (Int_t i=0; i<kSPECIES; i++) {
63     fTrackTime[i]=0;
64     fR[i]=0;
65     fITSr[i]=0;
66     fTPCr[i]=0;
67     fTRDr[i]=0;
68     fTOFr[i]=0;
69   }
70   Int_t i;
71   for (i=0; i<5; i++)   fRp[i]=0.;
72   for (i=0; i<15; i++)  fRc[i]=0.;
73   for (i=0; i<6; i++)   fITSindex[i]=0;
74   for (i=0; i<180; i++) fTPCindex[i]=0;
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     for (Int_t i=0;i<fTRDncls;i++) fTRDindex[i]=t->GetClusterIndex(i);
124     fTRDsignal=t->GetPIDsignal();
125     break;
126   default: 
127     Error("UpdateTrackParams()","Wrong flag !\n");
128     return kFALSE;
129   }
130
131   SetStatus(flags);
132   fLabel=t->GetLabel();
133
134   if (t->IsStartedTimeIntegral()) {
135     SetStatus(kTIME);
136     Double_t times[10];t->GetIntegratedTimes(times); SetIntegratedTimes(times);
137     SetIntegratedLength(t->GetIntegratedLength());
138   }
139
140   fRalpha=t->GetAlpha();
141   t->GetExternalParameters(fRx,fRp);
142   t->GetExternalCovariance(fRc);
143   
144   if (flags == kITSin)
145    {
146      AliKalmanTrack *itstrack = t;
147      if (itstrack)
148       {
149         itstrack->PropagateTo(3.,0.0028,65.19);
150         itstrack->PropagateToVertex();
151         
152         Double_t ralpha=t->GetAlpha();
153         Double_t rx;      // X-coordinate of the track reference plane 
154         Double_t rp[5];   // external track parameters  
155         t->GetExternalParameters(rx,rp);
156    
157         Double_t phi=TMath::ASin(rp[2]) + ralpha;
158         Double_t pt=1./TMath::Abs(rp[4]);
159         Double_t r=TMath::Sqrt(rx*rx + rp[0]*rp[0]);
160         
161         fVertexX=r*TMath::Cos(phi); 
162         fVertexY=r*TMath::Sin(phi); 
163         fVertexZ=rp[1]; 
164         
165         fVertexPx = pt*TMath::Cos(phi); 
166         fVertexPy = pt*TMath::Sin(phi); 
167         fVertexPz = pt*rp[3]; 
168         fVertex = kTRUE;
169       }
170    }
171   
172   return kTRUE;
173 }
174
175 //_______________________________________________________________________
176 void AliESDtrack::GetExternalParameters(Double_t &x, Double_t p[5]) const {
177   //---------------------------------------------------------------------
178   // This function returns external representation of the track parameters
179   //---------------------------------------------------------------------
180   x=fRx;
181   for (Int_t i=0; i<5; i++) p[i]=fRp[i];
182 }
183
184 Double_t AliESDtrack::GetP() const {
185   //---------------------------------------------------------------------
186   // This function returns the track momentum
187   //---------------------------------------------------------------------
188   Double_t lam=TMath::ATan(fRp[3]);
189   Double_t pt=1./TMath::Abs(fRp[4]);
190   return pt/TMath::Cos(lam);
191 }
192
193 void AliESDtrack::GetPxPyPz(Double_t *p) const {
194   //---------------------------------------------------------------------
195   // This function returns the global track momentum components
196   //---------------------------------------------------------------------
197   Double_t phi=TMath::ASin(fRp[2]) + fRalpha;
198   Double_t pt=1./TMath::Abs(fRp[4]);
199   p[0]=pt*TMath::Cos(phi); p[1]=pt*TMath::Sin(phi); p[2]=pt*fRp[3]; 
200 }
201
202 void AliESDtrack::GetXYZ(Double_t *xyz) const {
203   //---------------------------------------------------------------------
204   // This function returns the global track position
205   //---------------------------------------------------------------------
206   Double_t phi=TMath::ATan2(fRp[0],fRx) + fRalpha;
207   Double_t r=TMath::Sqrt(fRx*fRx + fRp[0]*fRp[0]);
208   xyz[0]=r*TMath::Cos(phi); xyz[1]=r*TMath::Sin(phi); xyz[2]=fRp[1]; 
209 }
210
211 //_______________________________________________________________________
212 void AliESDtrack::GetExternalCovariance(Double_t c[15]) const {
213   //---------------------------------------------------------------------
214   // This function returns external representation of the cov. matrix
215   //---------------------------------------------------------------------
216   for (Int_t i=0; i<15; i++) c[i]=fRc[i];
217 }
218
219 //_______________________________________________________________________
220 void AliESDtrack::GetIntegratedTimes(Double_t *times) const {
221   // Returns the array with integrated times for each particle hypothesis
222   for (Int_t i=0; i<kSPECIES; i++) times[i]=fTrackTime[i];
223 }
224
225 //_______________________________________________________________________
226 void AliESDtrack::SetIntegratedTimes(const Double_t *times) {
227   // Sets the array with integrated times for each particle hypotesis
228   for (Int_t i=0; i<kSPECIES; i++) fTrackTime[i]=times[i];
229 }
230
231 //_______________________________________________________________________
232 void AliESDtrack::SetITSpid(const Double_t *p) {
233   // Sets values for the probability of each particle type (in ITS)
234   for (Int_t i=0; i<kSPECIES; i++) fITSr[i]=p[i];
235   SetStatus(AliESDtrack::kITSpid);
236 }
237
238 //_______________________________________________________________________
239 void AliESDtrack::GetITSpid(Double_t *p) const {
240   // Gets the probability of each particle type (in ITS)
241   for (Int_t i=0; i<kSPECIES; i++) p[i]=fITSr[i];
242 }
243
244 //_______________________________________________________________________
245 Int_t AliESDtrack::GetITSclusters(UInt_t *idx) const {
246   //---------------------------------------------------------------------
247   // This function returns indices of the assgined ITS clusters 
248   //---------------------------------------------------------------------
249   for (Int_t i=0; i<fITSncls; i++) idx[i]=fITSindex[i];
250   return fITSncls;
251 }
252
253 //_______________________________________________________________________
254 Int_t AliESDtrack::GetTPCclusters(Int_t *idx) const {
255   //---------------------------------------------------------------------
256   // This function returns indices of the assgined ITS clusters 
257   //---------------------------------------------------------------------
258   for (Int_t i=0; i<180; i++) idx[i]=fTPCindex[i];  // MI I prefer some constant
259   return fTPCncls;
260 }
261
262 //_______________________________________________________________________
263 void AliESDtrack::SetTPCpid(const Double_t *p) {  
264   // Sets values for the probability of each particle type (in TPC)
265   for (Int_t i=0; i<kSPECIES; i++) fTPCr[i]=p[i];
266   SetStatus(AliESDtrack::kTPCpid);
267 }
268
269 //_______________________________________________________________________
270 void AliESDtrack::GetTPCpid(Double_t *p) const {
271   // Gets the probability of each particle type (in TPC)
272   for (Int_t i=0; i<kSPECIES; i++) p[i]=fTPCr[i];
273 }
274
275 //_______________________________________________________________________
276 Int_t AliESDtrack::GetTRDclusters(UInt_t *idx) const {
277   //---------------------------------------------------------------------
278   // This function returns indices of the assgined TRD clusters 
279   //---------------------------------------------------------------------
280   for (Int_t i=0; i<90; i++) idx[i]=fTRDindex[i];  // MI I prefer some constant
281   return fTRDncls;
282 }
283
284 //_______________________________________________________________________
285 void AliESDtrack::SetTRDpid(const Double_t *p) {  
286   // Sets values for the probability of each particle type (in TRD)
287   for (Int_t i=0; i<kSPECIES; i++) fTRDr[i]=p[i];
288   SetStatus(AliESDtrack::kTRDpid);
289 }
290
291 //_______________________________________________________________________
292 void AliESDtrack::GetTRDpid(Double_t *p) const {
293   // Gets the probability of each particle type (in TRD)
294   for (Int_t i=0; i<kSPECIES; i++) p[i]=fTRDr[i];
295 }
296
297 //_______________________________________________________________________
298 void    AliESDtrack::SetTRDpid(Int_t iSpecies, Float_t p)
299 {
300   // Sets the probability of particle type iSpecies to p (in TRD)
301   fTRDr[iSpecies] = p;
302 }
303
304 Float_t AliESDtrack::GetTRDpid(Int_t iSpecies) const
305 {
306   // Returns the probability of particle type iSpecies (in TRD)
307   return fTRDr[iSpecies];
308 }
309
310 //_______________________________________________________________________
311 void AliESDtrack::SetTOFpid(const Double_t *p) {  
312   // Sets the probability of each particle type (in TOF)
313   for (Int_t i=0; i<kSPECIES; i++) fTOFr[i]=p[i];
314   SetStatus(AliESDtrack::kTOFpid);
315 }
316
317 //_______________________________________________________________________
318 void AliESDtrack::GetTOFpid(Double_t *p) const {
319   // Gets probabilities of each particle type (in TOF)
320   for (Int_t i=0; i<kSPECIES; i++) p[i]=fTOFr[i];
321 }
322
323 //_______________________________________________________________________
324 void AliESDtrack::SetESDpid(const Double_t *p) {  
325   // Sets the probability of each particle type for the ESD track
326   for (Int_t i=0; i<kSPECIES; i++) fR[i]=p[i];
327   SetStatus(AliESDtrack::kESDpid);
328 }
329
330 //_______________________________________________________________________
331 void AliESDtrack::GetESDpid(Double_t *p) const {
332   // Gets probability of each particle type for the ESD track
333   for (Int_t i=0; i<kSPECIES; i++) p[i]=fR[i];
334 }
335
336 void AliESDtrack::GetVertexXYZ(Double_t& x,Double_t& y, Double_t&z) const
337 {
338 //returns track position in DCA to vertex  
339   x = fVertexX;
340   y = fVertexY;
341   z = fVertexZ;
342 }
343 void AliESDtrack::GetVertexPxPyPz(Double_t& px,Double_t& py, Double_t& pz) const
344 {
345 //returns track momentum in DCA to vertex  
346   px = fVertexPx;
347   py = fVertexPy;
348   pz = fVertexPz;
349 }