]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliESDtrack.h
Using symbolic constants instead of hardwired numbers (Yu.Belikov)
[u/mrichter/AliRoot.git] / STEER / AliESDtrack.h
1 #ifndef ALIESDTRACK_H
2 #define ALIESDTRACK_H
3 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4  * See cxx source for full Copyright notice                               */
5
6 /* $Id$ */
7
8 //-------------------------------------------------------------------------
9 //                          Class AliESDtrack
10 //   This is the class to deal with during the physics analysis of data
11 //      
12 //         Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch 
13 //-------------------------------------------------------------------------
14 /*****************************************************************************
15  *  Use GetExternalParameters() and GetExternalCovariance() to access the    *
16  *      track information regardless of its internal representation.         *
17  * This formation is now fixed in the following way:                         *
18  *      external param0:   local Y-coordinate of a track (cm)                *
19  *      external param1:   local Z-coordinate of a track (cm)                *
20  *      external param2:   local sine of the track momentum azimuthal angle  *
21  *      external param3:   tangent of the track momentum dip angle           *
22  *      external param4:   1/pt (1/(GeV/c))                                  *
23  *****************************************************************************/
24
25 #include <TBits.h>
26 #include "AliExternalTrackParam.h"
27 #include "AliPID.h"
28 #include <TVector3.h>
29
30 class AliESDVertex;
31 class AliKalmanTrack;
32 class AliTrackPointArray;
33
34 class AliESDtrack : public AliExternalTrackParam {
35 public:
36   AliESDtrack();
37   AliESDtrack(const AliESDtrack& track);
38   virtual ~AliESDtrack();
39   void MakeMiniESDtrack();
40   void SetID(Int_t id) { fID =id;}
41   Int_t GetID(){ return fID;}
42   void SetStatus(ULong_t flags) {fFlags|=flags;}
43   void ResetStatus(ULong_t flags) {fFlags&=~flags;}
44   Bool_t UpdateTrackParams(const AliKalmanTrack *t, ULong_t flags);
45   void SetIntegratedLength(Double_t l) {fTrackLength=l;}
46   void SetIntegratedTimes(const Double_t *times);
47   void SetESDpid(const Double_t *p);
48   void GetESDpid(Double_t *p) const;
49   
50   ULong_t GetStatus() const {return fFlags;}
51   Int_t GetLabel() const {return fLabel;}
52   void SetLabel(Int_t label) {fLabel = label;}
53
54   void GetExternalParameters(Double_t &x, Double_t p[5]) const;
55   void GetExternalCovariance(Double_t cov[15]) const;
56
57   Double_t GetIntegratedLength() const {return fTrackLength;}
58   void GetIntegratedTimes(Double_t *times) const;
59   Double_t GetMass() const;
60   TVector3 P3() const {Double_t p[3]; GetPxPyPz(p); return TVector3(p[0],p[1],p[2]);} //running track momentum
61   TVector3 X3() const {Double_t x[3]; GetXYZ(x); return TVector3(x[0],x[1],x[2]);}    //running track position 
62
63
64   Bool_t GetConstrainedPxPyPz(Double_t *p) const {
65     if (!fCp) return kFALSE;
66     return fCp->GetPxPyPz(p);
67   }
68   Bool_t GetConstrainedXYZ(Double_t *r) const {
69     if (!fCp) return kFALSE;
70     return fCp->GetXYZ(r);
71   }
72   Bool_t GetConstrainedExternalParameters
73               (Double_t &alpha, Double_t &x, Double_t p[5]) const;
74   Bool_t GetConstrainedExternalCovariance(Double_t cov[15]) const;
75   Double_t GetConstrainedChi2() const {return fCchi2;}
76
77
78   Bool_t GetInnerPxPyPz(Double_t *p) const {
79     if (!fIp) return kFALSE;
80     return fIp->GetPxPyPz(p);
81   }
82   const AliExternalTrackParam * GetInnerParam() const { return fIp;}
83   Bool_t GetInnerXYZ(Double_t *r) const {
84     if (!fIp) return kFALSE;
85     return fIp->GetXYZ(r);
86   }
87   Bool_t GetInnerExternalParameters
88         (Double_t &alpha, Double_t &x, Double_t p[5]) const;
89   Bool_t GetInnerExternalCovariance(Double_t cov[15]) const;
90  
91   const AliExternalTrackParam * GetOuterParam() const { return fOp;}
92   Bool_t GetOuterPxPyPz(Double_t *p) const {
93     if (!fOp) return kFALSE;
94     return fOp->GetPxPyPz(p);
95   }
96   Bool_t GetOuterXYZ(Double_t *r) const {
97     if (!fOp) return kFALSE;
98     return fOp->GetXYZ(r);
99   }
100   Bool_t GetOuterExternalParameters
101         (Double_t &alpha, Double_t &x, Double_t p[5]) const;
102   Bool_t GetOuterExternalCovariance(Double_t cov[15]) const;
103
104
105   Int_t GetNcls(Int_t idet) const;
106   Int_t GetClusters(Int_t idet, Int_t *idx) const;
107  
108   void SetITSpid(const Double_t *p);
109   void SetITSChi2MIP(const Float_t *chi2mip);
110   void SetITStrack(AliKalmanTrack * track){fITStrack=track;}
111   void GetITSpid(Double_t *p) const;
112   Float_t GetITSsignal() const {return fITSsignal;}
113   Float_t GetITSchi2() const {return fITSchi2;}
114   Int_t GetITSclusters(Int_t *idx) const;
115   Int_t GetITSLabel() const {return fITSLabel;}
116   Float_t GetITSFakeRatio() const {return fITSFakeRatio;}
117   AliKalmanTrack * GetITStrack(){return fITStrack;}
118
119   void SetTPCpid(const Double_t *p);
120   void GetTPCpid(Double_t *p) const;
121   void SetTPCPoints(Float_t points[4]){for (Int_t i=0;i<4;i++) fTPCPoints[i]=points[i];}
122   void SetTPCPointsF(UChar_t  findable){fTPCnclsF = findable;}
123   Float_t GetTPCPoints(Int_t i){return fTPCPoints[i];}
124   void SetKinkIndexes(Int_t points[3]) {for (Int_t i=0;i<3;i++) fKinkIndexes[i] = points[i];}
125   void SetV0Indexes(Int_t points[3]) {for (Int_t i=0;i<3;i++) fV0Indexes[i] = points[i];}
126   void SetTPCsignal(Float_t signal, Float_t sigma, UChar_t npoints){ fTPCsignal = signal; fTPCsignalS = sigma; fTPCsignalN = npoints;}
127   Float_t GetTPCsignal() const {return fTPCsignal;}
128   Float_t GetTPCchi2() const {return fTPCchi2;}
129   Int_t GetTPCclusters(Int_t *idx) const;
130   Float_t GetTPCdensity(Int_t row0, Int_t row1) const;
131   Int_t GetTPCLabel() const {return fTPCLabel;}
132   Int_t GetKinkIndex(Int_t i) const { return fKinkIndexes[i];}
133   Int_t GetV0Index(Int_t i) const { return fV0Indexes[i];}
134   const TBits& GetTPCClusterMap() const {return fTPCClusterMap;}
135   
136   void SetTRDpid(const Double_t *p);
137   void     SetTRDQuality(Float_t quality){fTRDQuality=quality;}
138   Float_t  GetTRDQuality()const {return fTRDQuality;}
139   void     SetTRDBudget(Float_t budget){fTRDBudget=budget;}
140   Float_t  GetTRDBudget()const {return fTRDBudget;}
141   void SetTRDtrack(AliKalmanTrack * track){fTRDtrack=track;}
142   void SetTRDsignals(Float_t dedx, Int_t i) {fTRDsignals[i]=dedx;}
143   void SetTRDTimBin(Int_t timbin, Int_t i) {fTRDTimBin[i]=timbin;}
144   void GetTRDpid(Double_t *p) const;
145   Float_t GetTRDsignal() const {return fTRDsignal;}
146   Float_t GetTRDsignals(Int_t i) const {return fTRDsignals[i];}
147   Int_t GetTRDTimBin(Int_t i) const {return fTRDTimBin[i];}
148   Float_t GetTRDchi2() const {return fTRDchi2;}
149   Int_t GetTRDclusters(Int_t *idx) const;
150   Int_t GetTRDncls() const {return fTRDncls;}
151   void    SetTRDpid(Int_t iSpecies, Float_t p);
152   Float_t GetTRDpid(Int_t iSpecies) const;
153   Int_t GetTRDLabel() const {return fTRDLabel;}
154
155
156   AliKalmanTrack * GetTRDtrack(){return fTRDtrack;}
157
158   void SetTOFsignal(Double_t tof) {fTOFsignal=tof;}
159   Float_t GetTOFsignal() const {return fTOFsignal;}
160   void SetTOFsignalToT(Double_t ToT) {fTOFsignalToT=ToT;}
161   Float_t GetTOFsignalToT() const {return fTOFsignalToT;}
162   Float_t GetTOFchi2() const {return fTOFchi2;}
163   void    SetTOFpid(const Double_t *p);
164   void    SetTOFLabel(const Int_t *p);
165   void    GetTOFpid(Double_t *p) const;
166   void    GetTOFLabel(Int_t *p) const;
167   void    GetTOFInfo(Float_t *info) const;
168   void    SetTOFInfo(Float_t *info);
169   Int_t   GetTOFCalChannel() const {return fTOFCalChannel;}
170   Int_t  GetTOFcluster() const {return fTOFindex;}
171   void  SetTOFcluster(Int_t index) {fTOFindex=index;}
172   void  SetTOFCalChannel(Int_t index) {fTOFCalChannel=index;}
173   
174   void    SetRICHsignal(Double_t beta) {fRICHsignal=beta;}
175   Float_t GetRICHsignal() const {return fRICHsignal;}
176   void    SetRICHpid(const Double_t *p);
177   void    GetRICHpid(Double_t *p) const;
178   void    SetRICHchi2(Double_t chi2) {fRICHchi2=chi2;}
179   Float_t GetRICHchi2() const {return fRICHchi2;}
180   void    SetRICHcluster(Int_t index) {fRICHindex=index;}
181   Int_t  GetRICHcluster() const {return fRICHindex;}
182   void    SetRICHnclusters(Int_t n) {fRICHncls=n;}
183   Int_t   GetRICHnclusters() const {return fRICHncls;}
184   void    SetRICHthetaPhi(Float_t  theta, Float_t  phi)      {fRICHtheta=theta; fRICHphi=phi;}
185   void    GetRICHthetaPhi(Float_t &theta, Float_t &phi)const {theta=fRICHtheta; phi=fRICHphi;}
186   void    SetRICHdxdy    (Float_t     dx, Float_t   dy)      {fRICHdx=dx;         fRICHdy=dy;}
187   void    GetRICHdxdy    (Float_t    &dx, Float_t  &dy)const {dx=fRICHdx;         dy=fRICHdy;}
188   void    SetRICHmipXY   (Float_t      x, Float_t    y)      {fRICHmipX=x;       fRICHmipY=y;} 
189   void    GetRICHmipXY   (Float_t     &x, Float_t   &y)const {x=fRICHmipX;       y=fRICHmipY;}
190   
191   Bool_t IsOn(Int_t mask) const {return (fFlags&mask)>0;}
192   Bool_t IsRICH()  const {return fFlags&kRICHpid;}
193   Bool_t IsPHOS()  const {return fFlags&kPHOSpid;}
194
195   void   SetTrackPointArray(AliTrackPointArray *points) { fPoints = points; }
196   AliTrackPointArray *GetTrackPointArray() const { return fPoints; }
197
198   Bool_t 
199     RelateToVertex(const AliESDVertex *vtx, Double_t b, Double_t maxd);
200   void GetImpactParameters(Float_t &xy,Float_t &z) const {xy=fD; z=fZ;}
201   void GetImpactParameters(Float_t p[2], Float_t cov[3]) const {
202     p[0]=fD; p[1]=fZ; cov[0]=fCdd; cov[1]=fCdz; cov[2]=fCzz;
203   }
204   virtual void Print(Option_t * opt) const ; 
205
206   enum {
207     kITSin=0x0001,kITSout=0x0002,kITSrefit=0x0004,kITSpid=0x0008,
208     kTPCin=0x0010,kTPCout=0x0020,kTPCrefit=0x0040,kTPCpid=0x0080,
209     kTRDin=0x0100,kTRDout=0x0200,kTRDrefit=0x0400,kTRDpid=0x0800,
210     kTOFin=0x1000,kTOFout=0x2000,kTOFrefit=0x4000,kTOFpid=0x8000,
211     kPHOSpid=0x10000, kRICHpid=0x20000,
212     kTRDbackup=0x80000,
213     kTRDStop=0x20000000,
214     kESDpid=0x40000000,
215     kTIME=0x80000000
216   }; 
217   enum {
218     kNPlane = 6,
219     kMaxITScluster=12,
220     kMaxTPCcluster=160,
221     kMaxTRDcluster=180
222   };
223 protected:
224   
225   //AliESDtrack & operator=(const AliESDtrack & );
226
227   ULong_t   fFlags;         // Reconstruction status flags 
228   Int_t     fLabel;         // Track label
229   Int_t     fID;            // Unique ID of the track
230   Float_t   fTrackLength;   // Track length
231   Float_t   fD;             // Impact parameter in XY plane
232   Float_t   fZ;             // Impact parameter in Z
233   Float_t   fCdd,fCdz,fCzz; // Covariance matrix of the impact parameters 
234   Float_t   fTrackTime[AliPID::kSPECIES]; // TOFs estimated by the tracking
235   Float_t   fR[AliPID::kSPECIES]; // combined "detector response probability"
236
237   Int_t   fStopVertex;  // Index of the stop vertex
238
239 //Track parameters constrained to the primary vertex
240   AliExternalTrackParam *fCp; 
241   Double_t fCchi2; //chi2 at the primary vertex
242
243 //Track parameters at the inner wall of the TPC
244   AliExternalTrackParam *fIp;
245
246 //Track parameters at the inner wall of the TRD 
247   AliExternalTrackParam *fOp;
248
249   // ITS related track information
250   Float_t fITSchi2;        // chi2 in the ITS
251   Float_t fITSchi2MIP[kMaxITScluster];     // chi2s in the ITS
252   Int_t   fITSncls;        // number of clusters assigned in the ITS
253   Int_t   fITSindex[kMaxITScluster];   //! indices of the assigned ITS clusters
254   Float_t fITSsignal;      // detector's PID signal
255   Float_t fITSr[AliPID::kSPECIES]; // "detector response probabilities" (for the PID)
256   Int_t   fITSLabel;       // label according TPC
257   Float_t fITSFakeRatio;   // ration of fake tracks
258   AliKalmanTrack * fITStrack; //! OWNER: pointer to the ITS track -- currently for debug purpose
259   
260   // TPC related track information
261   Float_t fTPCchi2;        // chi2 in the TPC
262   Int_t   fTPCncls;        // number of clusters assigned in the TPC
263   UShort_t fTPCnclsF;      // number of findable clusters in the TPC
264   Int_t  fTPCindex[kMaxTPCcluster];  //! indices of the assigned TPC clusters
265   TBits   fTPCClusterMap;  // Map of clusters, one bit per padrow; 1 if has a cluster on given padrow
266   Float_t fTPCsignal;      // detector's PID signal
267   UShort_t fTPCsignalN;      // number of points used for dEdx
268   Float_t  fTPCsignalS;    // RMS of dEdx measurement
269   Float_t fTPCr[AliPID::kSPECIES]; // "detector response probabilities" (for the PID)
270   Int_t   fTPCLabel;       // label according TPC
271   Float_t fTPCPoints[4];   // TPC points -first, max. dens, last and max density
272   Int_t   fKinkIndexes[3]; // array of indexes of posible kink candidates 
273   Int_t   fV0Indexes[3]; // array of indexes of posible kink candidates 
274
275   // TRD related track information
276   Float_t fTRDchi2;        // chi2 in the TRD
277   Int_t   fTRDncls;        // number of clusters assigned in the TRD
278   Int_t   fTRDncls0;       // number of clusters assigned in the TRD before first material cross
279   Int_t  fTRDindex[kMaxTRDcluster];   //! indices of the assigned TRD clusters
280   Float_t fTRDsignal;      // detector's PID signal
281   Float_t fTRDsignals[kNPlane];  // TRD signals from all six planes
282   Int_t   fTRDTimBin[kNPlane];   // Time bin of Max cluster from all six planes
283   Float_t fTRDr[AliPID::kSPECIES]; // "detector response probabilities" (for the PID)
284   Int_t   fTRDLabel;       // label according TRD
285   Float_t fTRDQuality;     //trd quality factor for TOF
286   Float_t fTRDBudget;     //trd material budget
287   AliKalmanTrack * fTRDtrack; //! OWNER: pointer to the TRD track -- currently for debug purpose
288
289   // TOF related track information
290   Float_t fTOFchi2;        // chi2 in the TOF
291   Int_t   fTOFindex;       // index of the assigned TOF cluster
292   Int_t   fTOFCalChannel; // Channel Index of the TOF Signal 
293   Float_t fTOFsignal;      // detector's PID signal
294   Float_t fTOFsignalToT;   // detector's ToT signal
295   Float_t fTOFr[AliPID::kSPECIES]; // "detector response probabilities" (for the PID)
296   Int_t   fTOFLabel[3];       // TOF label 
297   Float_t fTOFInfo[10];       //! TOF informations
298
299   // HMPID related track information
300   Float_t fRICHchi2;       // chi2 in the RICH
301   Int_t   fRICHncls;       // number of photon clusters
302   Int_t   fRICHindex;      // index of the assigned MIP cluster
303   Float_t fRICHsignal;     // RICH PID signal
304   Float_t fRICHr[AliPID::kSPECIES];// "detector response probabilities" (for the PID)
305   Float_t fRICHtheta;      // theta of the track extrapolated to the RICH
306   Float_t fRICHphi;        // phi of the track extrapolated to the RICH
307   Float_t fRICHdx;         // x of the track impact minus x of the MIP
308   Float_t fRICHdy;         // y of the track impact minus y of the MIP
309   Float_t fRICHmipX;       // x of the MIP in LORS
310   Float_t fRICHmipY;       // y of the MIP in LORS
311
312   AliTrackPointArray *fPoints; // Array which contains the track space points in the global frame
313
314   ClassDef(AliESDtrack,26)  //ESDtrack 
315 };
316
317 #endif 
318