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