]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliESDtrack.h
Possibility to set the track label
[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 <TObject.h>
27 #include "AliPID.h"
28 #include <TVector3.h>
29
30 class AliKalmanTrack;
31
32 const Int_t kNPlane = 6;
33
34 class AliESDtrack : public TObject {
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 SetImpactParameters(Float_t xy,Float_t z) {fD=xy; fZ=z;}
46   void SetIntegratedLength(Double_t l) {fTrackLength=l;}
47   void SetIntegratedTimes(const Double_t *times);
48   void SetESDpid(const Double_t *p);
49   void GetESDpid(Double_t *p) const;
50   
51   ULong_t GetStatus() const {return fFlags;}
52   Int_t GetLabel() const {return fLabel;}
53   void SetLabel(Int_t label) {fLabel = label;}
54   Double_t GetAlpha() const {return fRalpha;}
55   void GetExternalParameters(Double_t &x, Double_t p[5]) const;
56   void GetExternalCovariance(Double_t cov[15]) const;
57
58   Bool_t GetExternalParametersAt(Double_t x, Double_t p[5]) const;
59   Bool_t GetPxPyPzAt(Double_t x, Double_t p[3]) const;
60   Bool_t GetXYZAt(Double_t x, Double_t r[3]) const;
61
62   void GetImpactParameters(Float_t &xy,Float_t &z) const {xy=fD; z=fZ;}
63   Double_t GetIntegratedLength() const {return fTrackLength;}
64   void GetIntegratedTimes(Double_t *times) const;
65   Double_t GetMass() const;
66   Double_t GetP() const;
67   Bool_t GetPxPyPz(Double_t *p) const;
68   TVector3 P3() const {Double_t p[3]; GetPxPyPz(p); return TVector3(p[0],p[1],p[2]);} //running track momentum
69   Bool_t GetXYZ(Double_t *r) const;
70   TVector3 X3() const {Double_t x[3]; GetXYZ(x); return TVector3(x[0],x[1],x[2]);}    //running track position 
71   void GetCovariance(Double_t cov[21]) const;
72   Int_t GetSign() const {return (fRp[4]>0) ? 1 : -1;} 
73
74   void SetConstrainedTrackParams(const AliKalmanTrack *t, Double_t chi2);
75
76   Double_t GetConstrainedAlpha() const {return fCalpha;}
77   Double_t GetConstrainedChi2() const {return fCchi2;}
78   void GetConstrainedExternalParameters(Double_t &x, Double_t p[5]) const;
79   void GetConstrainedExternalCovariance(Double_t cov[15]) const;
80
81   Bool_t GetConstrainedPxPyPz(Double_t *p) const;
82   Bool_t GetConstrainedXYZ(Double_t *r) const;
83
84   Bool_t GetInnerPxPyPz(Double_t *p) const;
85   Bool_t GetInnerXYZ(Double_t *r) const;
86   void GetInnerExternalParameters(Double_t &x, Double_t p[5]) const;//skowron
87   void GetInnerExternalCovariance(Double_t cov[15]) const;//skowron
88   Double_t GetInnerAlpha() const {return fIalpha;}
89   
90   void SetITSpid(const Double_t *p);
91   void SetITSChi2MIP(const Float_t *chi2mip);
92   void SetITStrack(AliKalmanTrack * track){fITStrack=track;}
93   void GetITSpid(Double_t *p) const;
94   Float_t GetITSsignal() const {return fITSsignal;}
95   Float_t GetITSchi2() const {return fITSchi2;}
96   Int_t GetITSclusters(UInt_t *idx) const;
97   Int_t GetITSLabel() const {return fITSLabel;}
98   Float_t GetITSFakeRatio() const {return fITSFakeRatio;}
99   AliKalmanTrack * GetITStrack(){return fITStrack;}
100
101   void SetTPCpid(const Double_t *p);
102   void GetTPCpid(Double_t *p) const;
103   void SetTPCPoints(Float_t points[4]){for (Int_t i=0;i<4;i++) fTPCPoints[i]=points[i];}
104   void SetKinkIndexes(Int_t points[3]) {for (Int_t i=0;i<3;i++) fKinkIndexes[i] = points[i];}
105   void SetV0Indexes(Int_t points[3]) {for (Int_t i=0;i<3;i++) fV0Indexes[i] = points[i];}
106   Float_t GetTPCsignal() const {return fTPCsignal;}
107   Float_t GetTPCchi2() const {return fTPCchi2;}
108   Int_t GetTPCclusters(Int_t *idx) const;
109   Int_t GetTPCLabel() const {return fTPCLabel;}
110   Int_t GetKinkIndex(Int_t i) const { return fKinkIndexes[i];}
111   Int_t GetV0Index(Int_t i) const { return fV0Indexes[i];}
112   const TBits& GetTPCClusterMap() const {return fTPCClusterMap;}
113   
114   void SetTRDpid(const Double_t *p);
115   void     SetTRDQuality(Float_t quality){fTRDQuality=quality;}
116   Float_t  GetTRDQuality()const {return fTRDQuality;}
117   void SetTRDtrack(AliKalmanTrack * track){fTRDtrack=track;}
118   void SetTRDsignals(Float_t dedx, Int_t i) {fTRDsignals[i]=dedx;}
119   void SetTRDTimBin(Int_t timbin, Int_t i) {fTRDTimBin[i]=timbin;}
120   void GetTRDpid(Double_t *p) const;
121   Float_t GetTRDsignal() const {return fTRDsignal;}
122   Float_t GetTRDsignals(Int_t i) const {return fTRDsignals[i];}
123   Int_t GetTRDTimBin(Int_t i) const {return fTRDTimBin[i];}
124   Float_t GetTRDchi2() const {return fTRDchi2;}
125   Int_t GetTRDclusters(UInt_t *idx) const;
126   Int_t GetTRDncls() const {return fTRDncls;}
127   void    SetTRDpid(Int_t iSpecies, Float_t p);
128   Float_t GetTRDpid(Int_t iSpecies) const;
129   Int_t GetTRDLabel() const {return fTRDLabel;}
130   void GetTRDExternalParameters(Double_t &x, Double_t &alpha, Double_t p[5], Double_t cov[15]) const;//MI
131   AliKalmanTrack * GetTRDtrack(){return fTRDtrack;}
132
133   void SetTOFsignal(Double_t tof) {fTOFsignal=tof;}
134   Float_t GetTOFsignal() const {return fTOFsignal;}
135   Float_t GetTOFchi2() const {return fTOFchi2;}
136   void    SetTOFpid(const Double_t *p);
137   void    SetTOFLabel(const Int_t *p);
138   void    GetTOFpid(Double_t *p) const;
139   void    GetTOFLabel(Int_t *p) const;
140   void    GetTOFInfo(Float_t *info) const;
141   void    SetTOFInfo(Float_t *info);
142   UInt_t  GetTOFcluster() const {return fTOFindex;}
143   void  SetTOFcluster(UInt_t index) {fTOFindex=index;}
144   
145   void    SetRICHsignal(Double_t beta) {fRICHsignal=beta;}
146   Float_t GetRICHsignal() const {return fRICHsignal;}
147   void    SetRICHpid(const Double_t *p);
148   void    GetRICHpid(Double_t *p) const;
149   void    SetRICHchi2(Double_t chi2) {fRICHchi2=chi2;}
150   Float_t GetRICHchi2() const {return fRICHchi2;}
151   void    SetRICHcluster(UInt_t index) {fRICHindex=index;}
152   UInt_t  GetRICHcluster() const {return fRICHindex;}
153   void    SetRICHnclusters(Int_t n) {fRICHncls=n;}
154   Int_t   GetRICHnclusters() const {return fRICHncls;}
155   void    SetRICHthetaPhi(Double_t theta, Double_t phi) {
156     fRICHtheta=theta; fRICHphi=phi;
157   }
158   void    GetRICHthetaPhi(Double_t &theta, Double_t &phi) const {
159     theta=fRICHtheta; phi=fRICHphi;
160   }
161   void    SetRICHdxdy(Double_t dx, Double_t dy) {
162     fRICHdx=dx; fRICHdy=dy;
163   }
164   void    GetRICHdxdy(Double_t &dx, Double_t &dy) const {
165     dx=fRICHdx; dy=fRICHdy;
166   }
167   
168   void SetPHOSposition(const Double_t *pos)  {
169     fPHOSpos[0] = pos[0]; fPHOSpos[1]=pos[1]; fPHOSpos[2]=pos[2];
170   }
171   void SetPHOSsignal(Double_t ene) {fPHOSsignal = ene; }
172   void SetPHOSpid(const Double_t *p);
173   void GetPHOSposition(Double_t *pos) const {
174     pos[0]=fPHOSpos[0]; pos[1]=fPHOSpos[1]; pos[2]=fPHOSpos[2];
175   }
176   Float_t GetPHOSsignal() const {return fPHOSsignal;}
177   void GetPHOSpid(Double_t *p) const;  
178
179   void SetEMCALposition(const Double_t *pos)  {
180     fEMCALpos[0] = pos[0]; fEMCALpos[1]=pos[1]; fEMCALpos[2]=pos[2];
181   }
182   void SetEMCALsignal(Double_t ene) {fEMCALsignal = ene; }
183   void SetEMCALpid(const Double_t *p);
184   void GetEMCALposition(Double_t *pos) const {
185     pos[0]=fEMCALpos[0]; pos[1]=fEMCALpos[1]; pos[2]=fEMCALpos[2];
186   }
187   Float_t GetEMCALsignal() const {return fEMCALsignal;}
188   void GetEMCALpid(Double_t *p) const;  
189
190   Bool_t IsOn(Int_t mask) const {return (fFlags&mask)>0;}
191   Bool_t IsRICH()  const {return fFlags&kRICHpid;}
192   Bool_t IsPHOS()  const {return fFlags&kPHOSpid;}
193   Bool_t IsEMCAL() const {return fFlags&kEMCALpid;}
194
195   virtual void Print(Option_t * opt) const ; 
196
197   enum {
198     kITSin=0x0001,kITSout=0x0002,kITSrefit=0x0004,kITSpid=0x0008,
199     kTPCin=0x0010,kTPCout=0x0020,kTPCrefit=0x0040,kTPCpid=0x0080,
200     kTRDin=0x0100,kTRDout=0x0200,kTRDrefit=0x0400,kTRDpid=0x0800,
201     kTOFin=0x1000,kTOFout=0x2000,kTOFrefit=0x4000,kTOFpid=0x8000,
202     kPHOSpid=0x10000, kRICHpid=0x20000, kEMCALpid=0x40000,
203     kTRDbackup=0x80000,
204     kTRDStop=0x20000000,
205     kESDpid=0x40000000,
206     kTIME=0x80000000
207   }; 
208 protected:
209   
210   AliESDtrack & operator=(const AliESDtrack & );
211
212   ULong_t   fFlags;        // Reconstruction status flags 
213   Int_t     fLabel;        // Track label
214   Int_t     fID;           // Unique ID of the track
215   Float_t   fTrackLength;  // Track length
216   Float_t   fD;            // Impact parameter in XY-plane
217   Float_t   fZ;            // Impact parameter in Z 
218   Float_t   fTrackTime[AliPID::kSPECIES]; // TOFs estimated by the tracking
219   Float_t   fR[AliPID::kSPECIES];         // combined "detector response probability"
220
221   Int_t     fStopVertex;          // Index of stop vertex
222
223 //Running track parameters
224   Double_t fRalpha;  // track rotation angle
225   Double_t fRx;      // X-coordinate of the track reference plane 
226   Double_t fRp[5];   // external track parameters  
227   Double_t fRc[15];  // external cov. matrix of the track parameters
228
229 //Track parameters constrained to the primary vertex
230   Double_t fCalpha;   // Track rotation angle
231   Double_t fCx;       // x-coordinate of the track reference plane
232   Double_t fCp[5];    // external track parameters
233   Double_t fCc[15];   // external cov. matrix of the track parameters
234   Double_t fCchi2; //chi2 at the primary vertex
235
236 //Track parameters at the inner wall of the TPC
237   Double_t fIalpha;   // Track rotation angle
238   Double_t fIx;       // x-coordinate of the track reference plane
239   Double_t fIp[5];    // external track parameters
240   Double_t fIc[15];   // external cov. matrix of the track parameters
241
242 //Track parameters at the inner wall of the TRD 
243   Double_t fTalpha;   // Track rotation angle
244   Double_t fTx;       // x-coordinate of the track reference plane
245   Double_t fTp[5];    // external track parameters
246   Double_t fTc[15];   // external cov. matrix of the track parameters
247
248   // ITS related track information
249   Float_t fITSchi2;        // chi2 in the ITS
250   Float_t fITSchi2MIP[12];     // chi2s in the ITS
251   Int_t   fITSncls;        // number of clusters assigned in the ITS
252   UInt_t  fITSindex[6];    //! indices of the assigned ITS clusters
253   Float_t fITSsignal;      // detector's PID signal
254   Float_t fITSr[AliPID::kSPECIES]; // "detector response probabilities" (for the PID)
255   Int_t   fITSLabel;       // label according TPC
256   Float_t fITSFakeRatio;   // ration of fake tracks
257   AliKalmanTrack * fITStrack; //! OWNER: pointer to the ITS track -- currently for debug purpose
258   
259   // TPC related track information
260   Float_t fTPCchi2;        // chi2 in the TPC
261   Int_t   fTPCncls;        // number of clusters assigned in the TPC
262   Int_t  fTPCindex[180];  //! indices of the assigned TPC clusters
263   TBits   fTPCClusterMap;  // Map of clusters, one bit per padrow; 1 if has a cluster on given padrow
264   Float_t fTPCsignal;      // detector's PID signal
265   Float_t fTPCr[AliPID::kSPECIES]; // "detector response probabilities" (for the PID)
266   Int_t   fTPCLabel;       // label according TPC
267   Float_t fTPCPoints[4];   // TPC points -first, max. dens, last and max density
268   Int_t   fKinkIndexes[3]; // array of indexes of posible kink candidates 
269   Int_t   fV0Indexes[3]; // array of indexes of posible kink candidates 
270
271   // TRD related track information
272   Float_t fTRDchi2;        // chi2 in the TRD
273   Int_t   fTRDncls;        // number of clusters assigned in the TRD
274   Int_t   fTRDncls0;       // number of clusters assigned in the TRD before first material cross
275   UInt_t  fTRDindex[130];   //! indices of the assigned TRD clusters
276   Float_t fTRDsignal;      // detector's PID signal
277   Float_t fTRDsignals[kNPlane];  // TRD signals from all six planes
278   Int_t fTRDTimBin[kNPlane];     // Time bin of Max cluster from all six planes
279   Float_t fTRDr[AliPID::kSPECIES]; // "detector response probabilities" (for the PID)
280   Int_t   fTRDLabel;       // label according TRD
281   Float_t fTRDQuality;     //trd quality factor for TOF
282   AliKalmanTrack * fTRDtrack; //! OWNER: pointer to the TRD track -- currently for debug purpose
283
284   // TOF related track information
285   Float_t fTOFchi2;        // chi2 in the TOF
286   UInt_t  fTOFindex;       // index of the assigned TOF cluster
287   Float_t fTOFsignal;      // detector's PID signal
288   Float_t fTOFr[AliPID::kSPECIES]; // "detector response probabilities" (for the PID)
289   Int_t   fTOFLabel[3];       // TOF label 
290   Float_t fTOFInfo[10];       //! TOF informations
291
292   // PHOS related track information 
293   Float_t fPHOSpos[3]; // position localised by PHOS in global coordinate system
294   Float_t fPHOSsignal; // energy measured by PHOS
295   Float_t fPHOSr[AliPID::kSPECIESN]; // PID information from PHOS
296
297   // EMCAL related track information 
298   Float_t fEMCALpos[3]; //position localised by EMCAL in global coordinate system
299   Float_t fEMCALsignal; // energy measured by EMCAL
300   Float_t fEMCALr[AliPID::kSPECIESN]; // PID information from EMCAL
301
302   // HMPID related track information
303   Float_t fRICHchi2;       // chi2 in the RICH
304   Int_t   fRICHncls;       // number of photon clusters
305   UInt_t  fRICHindex;      // index of the assigned MIP cluster
306   Float_t fRICHsignal;     // RICH PID signal
307   Float_t fRICHr[AliPID::kSPECIES];// "detector response probabilities" (for the PID)
308   Float_t fRICHtheta;      // theta of the track extrapolated to the RICH
309   Float_t fRICHphi;        // phi of the track extrapolated to the RICH
310   Float_t fRICHdx;         // x of the track impact minus x of the MIP
311   Float_t fRICHdy;         // y of the track impact minus y of the MIP
312         
313   ClassDef(AliESDtrack,13)  //ESDtrack 
314 };
315
316 #endif 
317