e6326636fad73b95c742269bd46ce4fc3f8b7c2d
[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      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::GetExternalParametersAt(Double_t x, Double_t p[5]) const {
177   //---------------------------------------------------------------------
178   // This function returns external representation of the track parameters
179   // at the plane x
180   //---------------------------------------------------------------------
181   Double_t dx=x-fRx;
182   Double_t c=fRp[4]/AliKalmanTrack::GetConvConst();
183   Double_t f1=fRp[2], f2=f1 + c*dx;
184   Double_t r1=sqrt(1.- f1*f1), r2=sqrt(1.- f2*f2);    
185
186   p[0]=fRp[0]+dx*(f1+f2)/(r1+r2);
187   p[1]=fRp[1]+dx*(f1+f2)/(f1*r2 + f2*r1)*fRp[3];
188   p[2]=fRp[2]+dx*c;
189   p[3]=fRp[3];
190   p[4]=fRp[4];
191 }
192
193 //_______________________________________________________________________
194 void AliESDtrack::GetExternalParameters(Double_t &x, Double_t p[5]) const {
195   //---------------------------------------------------------------------
196   // This function returns external representation of the track parameters
197   //---------------------------------------------------------------------
198   x=fRx;
199   for (Int_t i=0; i<5; i++) p[i]=fRp[i];
200 }
201
202 Double_t AliESDtrack::GetP() const {
203   //---------------------------------------------------------------------
204   // This function returns the track momentum
205   //---------------------------------------------------------------------
206   Double_t lam=TMath::ATan(fRp[3]);
207   Double_t pt=1./TMath::Abs(fRp[4]);
208   return pt/TMath::Cos(lam);
209 }
210
211 void AliESDtrack::GetPxPyPz(Double_t *p) const {
212   //---------------------------------------------------------------------
213   // This function returns the global track momentum components
214   //---------------------------------------------------------------------
215   Double_t phi=TMath::ASin(fRp[2]) + fRalpha;
216   Double_t pt=1./TMath::Abs(fRp[4]);
217   p[0]=pt*TMath::Cos(phi); p[1]=pt*TMath::Sin(phi); p[2]=pt*fRp[3]; 
218 }
219
220 void AliESDtrack::GetXYZ(Double_t *xyz) const {
221   //---------------------------------------------------------------------
222   // This function returns the global track position
223   //---------------------------------------------------------------------
224   Double_t phi=TMath::ASin(fRp[2]) + fRalpha;
225   Double_t r=TMath::Sqrt(fRx*fRx + fRp[0]*fRp[0]);
226   xyz[0]=r*TMath::Cos(phi); xyz[1]=r*TMath::Sin(phi); xyz[2]=fRp[1]; 
227 }
228
229 //_______________________________________________________________________
230 void AliESDtrack::GetExternalCovariance(Double_t c[15]) const {
231   //---------------------------------------------------------------------
232   // This function returns external representation of the cov. matrix
233   //---------------------------------------------------------------------
234   for (Int_t i=0; i<15; i++) c[i]=fRc[i];
235 }
236
237 //_______________________________________________________________________
238 void AliESDtrack::GetIntegratedTimes(Double_t *times) const {
239   // Returns the array with integrated times for each particle hypothesis
240   for (Int_t i=0; i<kSPECIES; i++) times[i]=fTrackTime[i];
241 }
242
243 //_______________________________________________________________________
244 void AliESDtrack::SetIntegratedTimes(const Double_t *times) {
245   // Sets the array with integrated times for each particle hypotesis
246   for (Int_t i=0; i<kSPECIES; i++) fTrackTime[i]=times[i];
247 }
248
249 //_______________________________________________________________________
250 void AliESDtrack::SetITSpid(const Double_t *p) {
251   // Sets values for the probability of each particle type (in ITS)
252   for (Int_t i=0; i<kSPECIES; i++) fITSr[i]=p[i];
253   SetStatus(AliESDtrack::kITSpid);
254 }
255
256 //_______________________________________________________________________
257 void AliESDtrack::GetITSpid(Double_t *p) const {
258   // Gets the probability of each particle type (in ITS)
259   for (Int_t i=0; i<kSPECIES; i++) p[i]=fITSr[i];
260 }
261
262 //_______________________________________________________________________
263 Int_t AliESDtrack::GetITSclusters(UInt_t *idx) const {
264   //---------------------------------------------------------------------
265   // This function returns indices of the assgined ITS clusters 
266   //---------------------------------------------------------------------
267   for (Int_t i=0; i<fITSncls; i++) idx[i]=fITSindex[i];
268   return fITSncls;
269 }
270
271 //_______________________________________________________________________
272 Int_t AliESDtrack::GetTPCclusters(Int_t *idx) const {
273   //---------------------------------------------------------------------
274   // This function returns indices of the assgined ITS clusters 
275   //---------------------------------------------------------------------
276   for (Int_t i=0; i<180; i++) idx[i]=fTPCindex[i];  // MI I prefer some constant
277   return fTPCncls;
278 }
279
280 //_______________________________________________________________________
281 void AliESDtrack::SetTPCpid(const Double_t *p) {  
282   // Sets values for the probability of each particle type (in TPC)
283   for (Int_t i=0; i<kSPECIES; i++) fTPCr[i]=p[i];
284   SetStatus(AliESDtrack::kTPCpid);
285 }
286
287 //_______________________________________________________________________
288 void AliESDtrack::GetTPCpid(Double_t *p) const {
289   // Gets the probability of each particle type (in TPC)
290   for (Int_t i=0; i<kSPECIES; i++) p[i]=fTPCr[i];
291 }
292
293 //_______________________________________________________________________
294 void AliESDtrack::SetTRDpid(const Double_t *p) {  
295   // Sets values for the probability of each particle type (in TRD)
296   for (Int_t i=0; i<kSPECIES; i++) fTRDr[i]=p[i];
297   SetStatus(AliESDtrack::kTRDpid);
298 }
299
300 //_______________________________________________________________________
301 void AliESDtrack::GetTRDpid(Double_t *p) const {
302   // Gets the probability of each particle type (in TRD)
303   for (Int_t i=0; i<kSPECIES; i++) p[i]=fTRDr[i];
304 }
305
306 //_______________________________________________________________________
307 void    AliESDtrack::SetTRDpid(Int_t iSpecies, Float_t p)
308 {
309   // Sets the probability of particle type iSpecies to p (in TRD)
310   fTRDr[iSpecies] = p;
311 }
312
313 Float_t AliESDtrack::GetTRDpid(Int_t iSpecies) const
314 {
315   // Returns the probability of particle type iSpecies (in TRD)
316   return fTRDr[iSpecies];
317 }
318
319 //_______________________________________________________________________
320 void AliESDtrack::SetTOFpid(const Double_t *p) {  
321   // Sets the probability of each particle type (in TOF)
322   for (Int_t i=0; i<kSPECIES; i++) fTOFr[i]=p[i];
323   SetStatus(AliESDtrack::kTOFpid);
324 }
325
326 //_______________________________________________________________________
327 void AliESDtrack::GetTOFpid(Double_t *p) const {
328   // Gets probabilities of each particle type (in TOF)
329   for (Int_t i=0; i<kSPECIES; i++) p[i]=fTOFr[i];
330 }
331
332 //_______________________________________________________________________
333 void AliESDtrack::SetESDpid(const Double_t *p) {  
334   // Sets the probability of each particle type for the ESD track
335   for (Int_t i=0; i<kSPECIES; i++) fR[i]=p[i];
336   SetStatus(AliESDtrack::kESDpid);
337 }
338
339 //_______________________________________________________________________
340 void AliESDtrack::GetESDpid(Double_t *p) const {
341   // Gets probability of each particle type for the ESD track
342   for (Int_t i=0; i<kSPECIES; i++) p[i]=fR[i];
343 }
344
345 void AliESDtrack::GetVertexXYZ(Double_t& x,Double_t& y, Double_t&z) const
346 {
347 //returns track position in DCA to vertex  
348   x = fVertexX;
349   y = fVertexY;
350   z = fVertexZ;
351 }
352 void AliESDtrack::GetVertexPxPyPz(Double_t& px,Double_t& py, Double_t& pz) const
353 {
354 //returns track momentum in DCA to vertex  
355   px = fVertexPx;
356   py = fVertexPy;
357   pz = fVertexPz;
358 }