]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliESDtrack.h
Fixes for Coverity warnings - M. Sitta
[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  * The Get*Label() getters return the label of the associated MC particle.   *
25  * The absolute value of this label is the index of the particle within the  *
26  * MC stack. If the label is negative, this track was assigned a certain     *
27  * number of clusters that did not in fact belong to this track.             *
28  *****************************************************************************/
29
30 #include <TBits.h>
31 #include "AliExternalTrackParam.h"
32 #include "AliVTrack.h"
33 #include "AliPID.h"
34 #include "AliESDfriendTrack.h"
35
36 class TParticle;
37 class AliESDVertex;
38 class AliESDEvent;
39 class AliKalmanTrack;
40 class AliTrackPointArray;
41 class TPolyMarker3D;
42
43 class AliESDtrack : public AliExternalTrackParam {
44 public:
45   enum {
46     kITSin=0x0001,kITSout=0x0002,kITSrefit=0x0004,kITSpid=0x0008,
47     kTPCin=0x0010,kTPCout=0x0020,kTPCrefit=0x0040,kTPCpid=0x0080,
48     kTRDin=0x0100,kTRDout=0x0200,kTRDrefit=0x0400,kTRDpid=0x0800,
49     kTOFin=0x1000,kTOFout=0x2000,kTOFrefit=0x4000,kTOFpid=0x8000,
50     kTOFmismatch=0x100000,
51     kHMPIDout=0x10000,kHMPIDpid=0x20000,
52     kEMCALmatch=0x40000,
53     kPHOSmatch=0x200000,
54     kTRDbackup =0x80000,
55     kTRDStop=0x20000000,
56     kESDpid=0x40000000,
57     kTIME=0x80000000,
58     kGlobalMerge=0x08000000,
59     kITSpureSA=0x10000000,
60     kMultInV0=0x2000000,    //BIT(25): assumed to be belong to V0 in multiplicity estimates
61     kMultSec=0x4000000     //BIT(26): assumed to be secondary (due to the DCA) in multiplicity estimates
62   }; 
63   enum {
64     kTRDnPlanes = 6,
65     kEMCALNoMatch = -4096
66   };
67   AliESDtrack();
68   AliESDtrack(const AliESDtrack& track);
69   AliESDtrack(const AliVTrack* track);
70   AliESDtrack(TParticle * part);
71   virtual ~AliESDtrack();
72   virtual void Copy(TObject &obj) const;
73   const AliESDfriendTrack *GetFriendTrack() const {return fFriendTrack;}
74   void SetFriendTrack(const AliESDfriendTrack *t) {
75     delete fFriendTrack; fFriendTrack=new AliESDfriendTrack(*t);
76   }
77   void ReleaseESDfriendTrack() { delete fFriendTrack;  fFriendTrack=0; }
78   void AddCalibObject(TObject * object);     // add calib object to the list
79   TObject *  GetCalibObject(Int_t index);    // return calib objct at given position
80   void MakeMiniESDtrack();
81   void SetID(Short_t id) { fID =id;}
82   Int_t GetID() const { return fID;}
83   void SetVertexID(Char_t id) { fVertexID=id;}
84   Char_t GetVertexID() const { return fVertexID;}
85   void SetStatus(ULong_t flags) {fFlags|=flags;}
86   void ResetStatus(ULong_t flags) {fFlags&=~flags;}
87   Bool_t UpdateTrackParams(const AliKalmanTrack *t, ULong_t flags);
88   void SetIntegratedLength(Double_t l) {fTrackLength=l;}
89   void SetIntegratedTimes(const Double_t *times);
90   void SetESDpid(const Double_t *p);
91   void GetESDpid(Double_t *p) const;
92   virtual const Double_t *PID() const { return fR; }
93
94   Bool_t IsOn(Int_t mask) const {return (fFlags&mask)>0;}
95   ULong_t GetStatus() const {return fFlags;}
96   Int_t GetLabel() const {return fLabel;}
97   void SetLabel(Int_t label) {fLabel = label;}
98
99   void GetExternalParameters(Double_t &x, Double_t p[5]) const;
100   void GetExternalCovariance(Double_t cov[15]) const;
101
102   Double_t GetIntegratedLength() const {return fTrackLength;}
103   void GetIntegratedTimes(Double_t *times) const;
104   Double_t GetMass() const;
105   Double_t M() const;
106   Double_t E() const;
107   Double_t Y() const;
108
109   Bool_t GetConstrainedPxPyPz(Double_t *p) const {
110     if (!fCp) return kFALSE;
111     return fCp->GetPxPyPz(p);
112   }
113   Bool_t GetConstrainedXYZ(Double_t *r) const {
114     if (!fCp) return kFALSE;
115     return fCp->GetXYZ(r);
116   }
117   const AliExternalTrackParam *GetConstrainedParam() const {return fCp;}
118   Bool_t GetConstrainedExternalParameters
119               (Double_t &alpha, Double_t &x, Double_t p[5]) const;
120   Bool_t GetConstrainedExternalCovariance(Double_t cov[15]) const;
121   Double_t GetConstrainedChi2() const {return fCchi2;}
122   //
123   
124   // global track chi2
125   void SetGlobalChi2(Double_t chi2) {fGlobalChi2 = chi2;}
126   Double_t GetGlobalChi2() const {return fGlobalChi2;}
127
128   Bool_t GetInnerPxPyPz(Double_t *p) const {
129     if (!fIp) return kFALSE;
130     return fIp->GetPxPyPz(p);
131   }
132   const AliExternalTrackParam * GetInnerParam() const { return fIp;}
133   const AliExternalTrackParam * GetTPCInnerParam() const {return fTPCInner;}
134   Bool_t FillTPCOnlyTrack(AliESDtrack &track);
135   Bool_t GetInnerXYZ(Double_t *r) const {
136     if (!fIp) return kFALSE;
137     return fIp->GetXYZ(r);
138   }
139   Bool_t GetInnerExternalParameters
140         (Double_t &alpha, Double_t &x, Double_t p[5]) const;
141   Bool_t GetInnerExternalCovariance(Double_t cov[15]) const;
142  
143   void SetOuterParam(const AliExternalTrackParam *p, ULong_t flags);
144
145   void SetOuterHmpParam(const AliExternalTrackParam *p, ULong_t flags);
146
147   const AliExternalTrackParam * GetOuterParam() const { return fOp;}
148
149   const AliExternalTrackParam * GetOuterHmpParam() const { return fHMPIDp;}
150   
151   Bool_t GetOuterPxPyPz(Double_t *p) const {
152     if (!fOp) return kFALSE;
153     return fOp->GetPxPyPz(p);
154   }
155   Bool_t GetOuterHmpPxPyPz(Double_t *p) const {
156     if (!fHMPIDp) return kFALSE;
157     return fHMPIDp->GetPxPyPz(p);
158   }
159   
160   Bool_t GetOuterXYZ(Double_t *r) const {
161     if (!fOp) return kFALSE;
162     return fOp->GetXYZ(r);
163   }
164     Bool_t GetOuterHmpXYZ(Double_t *r) const {
165     if (!fHMPIDp) return kFALSE;
166     return fHMPIDp->GetXYZ(r);
167   }
168
169   Bool_t GetOuterExternalParameters
170         (Double_t &alpha, Double_t &x, Double_t p[5]) const;
171   Bool_t GetOuterExternalCovariance(Double_t cov[15]) const;
172
173   Bool_t GetOuterHmpExternalParameters
174         (Double_t &alpha, Double_t &x, Double_t p[5]) const;
175   Bool_t GetOuterHmpExternalCovariance(Double_t cov[15]) const;
176
177   
178   Int_t GetNcls(Int_t idet) const;
179   Int_t GetClusters(Int_t idet, Int_t *idx) const;
180  
181   void    SetITSpid(const Double_t *p);
182   void    GetITSpid(Double_t *p) const;
183
184   Double_t GetITSsignal() const {return fITSsignal;}
185   void    SetITSdEdxSamples(const Double_t s[4]);
186   void    GetITSdEdxSamples(Double_t *s) const;
187
188   Double_t GetITSchi2() const {return fITSchi2;}
189   Char_t   GetITSclusters(Int_t *idx) const;
190   UChar_t GetITSClusterMap() const {return fITSClusterMap;}
191   void    SetITSModuleIndex(Int_t ilayer,Int_t idx) {fITSModule[ilayer]=idx;}
192   Int_t   GetITSModuleIndex(Int_t ilayer) const {return fITSModule[ilayer];}
193   Bool_t  GetITSModuleIndexInfo(Int_t ilayer,Int_t &idet,Int_t &status,
194                                 Float_t &xloc,Float_t &zloc) const;
195   Int_t   GetITSLabel() const {return fITSLabel;}
196   void    SetITStrack(AliKalmanTrack * track){
197      fFriendTrack->SetITStrack(track);
198   }
199   AliKalmanTrack *GetITStrack(){
200      return fFriendTrack->GetITStrack();
201   }
202   Bool_t  HasPointOnITSLayer(Int_t i) const {return TESTBIT(fITSClusterMap,i);}
203
204   void    SetTPCpid(const Double_t *p);
205   void    GetTPCpid(Double_t *p) const;
206   void    SetTPCPoints(Float_t points[4]){
207      for (Int_t i=0;i<4;i++) fTPCPoints[i]=points[i];
208   }
209   void    SetTPCPointsF(UChar_t  findable){fTPCnclsF = findable;}
210   void    SetTPCPointsFIter1(UChar_t  findable){fTPCnclsFIter1 = findable;}
211   UShort_t   GetTPCNcls() const { return fTPCncls;}
212   UShort_t   GetTPCNclsF() const { return fTPCnclsF;}
213   UShort_t   GetTPCNclsIter1() const { return fTPCnclsIter1;}
214   UShort_t   GetTPCNclsFIter1() const { return fTPCnclsFIter1;}
215   UShort_t   GetTPCnclsS(Int_t i0=0,Int_t i1=159) const;
216   UShort_t   GetTPCncls(Int_t row0=0,Int_t row1=159) const;
217   Double_t GetTPCPoints(Int_t i) const {return fTPCPoints[i];}
218   void    SetKinkIndexes(Int_t points[3]) {
219      for (Int_t i=0;i<3;i++) fKinkIndexes[i] = points[i];
220   }
221   void    SetV0Indexes(Int_t points[3]) {
222      for (Int_t i=0;i<3;i++) fV0Indexes[i] = points[i];
223   }
224   void    SetTPCsignal(Float_t signal, Float_t sigma, UChar_t npoints){ 
225      fTPCsignal = signal; fTPCsignalS = sigma; fTPCsignalN = npoints;
226   }
227   Double_t GetTPCsignal() const {return fTPCsignal;}
228   Double_t GetTPCsignalSigma() const {return fTPCsignalS;}
229   UShort_t GetTPCsignalN() const {return fTPCsignalN;}
230   Double_t GetTPCchi2() const {return fTPCchi2;}
231   Double_t GetTPCchi2Iter1() const {return fTPCchi2Iter1;}
232   UShort_t   GetTPCclusters(Int_t *idx) const;
233   Double_t GetTPCdensity(Int_t row0, Int_t row1) const;
234   Int_t   GetTPCLabel() const {return fTPCLabel;}
235   Int_t   GetKinkIndex(Int_t i) const { return fKinkIndexes[i];}
236   Int_t   GetV0Index(Int_t i) const { return fV0Indexes[i];}
237   const TBits& GetTPCClusterMap() const {return fTPCClusterMap;}
238   const TBits& GetTPCSharedMap() const {return fTPCSharedMap;}
239   void    SetTPCClusterMap(const TBits amap) {fTPCClusterMap = amap;}
240   void    SetTPCSharedMap(const TBits amap) {fTPCSharedMap = amap;}
241   Float_t GetTPCClusterInfo(Int_t nNeighbours=3, Int_t type=0) const;
242   
243   void    SetTRDpid(const Double_t *p);
244   
245 // A.Bercuci
246   void    SetTRDntracklets(UChar_t q){fTRDntracklets = q;}
247   UChar_t GetTRDntracklets() const {return (fTRDntracklets>>3)&7;}
248   UChar_t GetTRDntrackletsPID() const {return fTRDntracklets&7;}
249   // TEMPORARY alias asked by the HFE group to allow 
250   // reading of the v4-16-Release data with TRUNK related software (A.Bercuci@Apr 30th 09) 
251   UChar_t GetTRDpidQuality() const {return GetTRDntrackletsPID();}
252 // end A.Bercuci
253
254   void     SetNumberOfTRDslices(Int_t n);
255   Int_t    GetNumberOfTRDslices() const;
256   void     SetTRDslice(Double_t q, Int_t plane, Int_t slice);
257   void     SetTRDmomentum(Double_t p, Int_t plane, Double_t *sp=0x0);
258   Double_t GetTRDslice(Int_t plane, Int_t slice=-1) const;
259   Double_t GetTRDmomentum(Int_t plane, Double_t *sp=0x0) const;
260         
261   void    SetTRDQuality(Float_t quality){fTRDQuality=quality;}
262   Double_t GetTRDQuality()const {return fTRDQuality;}
263   void    SetTRDBudget(Float_t budget){fTRDBudget=budget;}
264   Double_t GetTRDBudget()const {return fTRDBudget;}
265
266   void    SetTRDTimBin(Int_t timbin, Int_t i) {fTRDTimBin[i]=timbin;}
267   void    GetTRDpid(Double_t *p) const;
268   Double_t GetTRDsignal() const {return fTRDsignal;}
269
270   Char_t   GetTRDTimBin(Int_t i) const {return fTRDTimBin[i];}
271   Double_t GetTRDchi2() const {return fTRDchi2;}
272   UChar_t   GetTRDclusters(Int_t *idx) const;
273   UChar_t   GetTRDncls() const {return fTRDncls;}
274   UChar_t   GetTRDncls0() const {return fTRDncls0;}
275   UChar_t   GetTRDtracklets(Int_t *idx) const;
276   void    SetTRDpid(Int_t iSpecies, Float_t p);
277   Double_t GetTRDpid(Int_t iSpecies) const;
278   Int_t   GetTRDLabel() const {return fTRDLabel;}
279
280   void    SetTRDtrack(AliKalmanTrack * track){
281      fFriendTrack->SetTRDtrack(track);
282   }
283   AliKalmanTrack *GetTRDtrack(){
284      return fFriendTrack->GetTRDtrack();
285   }
286
287   void    SetTOFsignal(Double_t tof) {fTOFsignal=tof;}
288   Double_t GetTOFsignal() const {return fTOFsignal;}
289   void    SetTOFsignalToT(Double_t ToT) {fTOFsignalToT=ToT;}
290   Double_t GetTOFsignalToT() const {return fTOFsignalToT;}
291   void    SetTOFsignalRaw(Double_t tof) {fTOFsignalRaw=tof;}
292   Double_t GetTOFsignalRaw() const {return fTOFsignalRaw;}
293   void    SetTOFsignalDz(Double_t dz) {fTOFsignalDz=dz;}
294   Double_t GetTOFsignalDz() const {return fTOFsignalDz;}
295   void    SetTOFsignalDx(Double_t dx) {fTOFsignalDx=dx;}
296   Double_t GetTOFsignalDx() const {return fTOFsignalDx;}
297   void     SetTOFDeltaBC(Short_t deltaBC) {fTOFdeltaBC=deltaBC;};
298   Short_t  GetTOFDeltaBC() const {return fTOFdeltaBC;}
299   void     SetTOFL0L1(Short_t l0l1) {fTOFl0l1=l0l1;};
300   Short_t  GetTOFL0L1() const {return fTOFl0l1;}
301   Double_t GetTOFchi2() const {return fTOFchi2;}
302   void    SetTOFpid(const Double_t *p);
303   void    SetTOFLabel(const Int_t *p);
304   void    GetTOFpid(Double_t *p) const;
305   void    GetTOFLabel(Int_t *p) const;
306   void    GetTOFInfo(Float_t *info) const;
307   void    SetTOFInfo(Float_t *info);
308   Int_t   GetTOFCalChannel() const {return fTOFCalChannel;}
309   Int_t   GetTOFcluster() const {return fTOFindex;}
310   void    SetTOFcluster(Int_t index) {fTOFindex=index;}
311   void    SetTOFCalChannel(Int_t index) {fTOFCalChannel=index;}
312
313 // HMPID methodes +++++++++++++++++++++++++++++++++ (kir)
314   void    SetHMPIDsignal(Double_t theta) {fHMPIDsignal=theta;}
315   Double_t GetHMPIDsignal() const {return fHMPIDsignal;}
316   void    SetHMPIDpid(const Double_t *p);
317   void    GetHMPIDpid(Double_t *p) const;  
318   void    SetHMPIDchi2(Double_t chi2) {fHMPIDchi2=chi2;}
319   Double_t GetHMPIDchi2() const {return fHMPIDchi2;}
320   void    SetHMPIDcluIdx(Int_t ch,Int_t idx) {fHMPIDcluIdx=ch*1000000+idx;}
321   Int_t   GetHMPIDcluIdx() const {return fHMPIDcluIdx;}
322   void    SetHMPIDtrk(Float_t  x, Float_t  y, Float_t  th, Float_t  ph) {
323      fHMPIDtrkX=x; fHMPIDtrkY=y; fHMPIDtrkTheta=th; fHMPIDtrkPhi=ph;
324   }
325   void    GetHMPIDtrk(Float_t &x, Float_t &y, Float_t &th, Float_t &ph) const {
326      x=fHMPIDtrkX; y=fHMPIDtrkY; th=fHMPIDtrkTheta; ph=fHMPIDtrkPhi;
327   }
328   void    SetHMPIDmip(Float_t  x, Float_t  y, Int_t q, Int_t nph=0) {
329      fHMPIDmipX=x; fHMPIDmipY=y; fHMPIDqn=1000000*nph+q;
330   }
331   void    GetHMPIDmip(Float_t &x,Float_t &y,Int_t &q,Int_t &nph) const {
332      x=fHMPIDmipX; y=fHMPIDmipY; q=fHMPIDqn%1000000; nph=fHMPIDqn/1000000;
333   }
334   Bool_t  IsHMPID() const {return fFlags&kHMPIDpid;}
335   Bool_t  IsPureITSStandalone() const {return fFlags&kITSpureSA;}
336   Bool_t  IsMultPrimary() const {return !(fFlags&kMultSec);}
337   Bool_t  IsMultSecondary() const {return (fFlags&kMultSec);}
338
339   Int_t GetEMCALcluster() {return fCaloIndex;}
340   void SetEMCALcluster(Int_t index) {fCaloIndex=index;}
341   Bool_t IsEMCAL() const {return fFlags&kEMCALmatch;}
342
343   Int_t GetPHOScluster() {return fCaloIndex;}
344   void SetPHOScluster(Int_t index) {fCaloIndex=index;}
345   Bool_t IsPHOS() const {return fFlags&kPHOSmatch;}
346   Double_t GetPHOSdx()const{return fCaloDx ;}
347   Double_t GetPHOSdz()const{return fCaloDz ;}
348   void SetPHOSdxdz(Double_t dx, Double_t dz){fCaloDx=dx,fCaloDz=dz;}
349
350
351   void SetTrackPointArray(AliTrackPointArray *points) {
352     fFriendTrack->SetTrackPointArray(points);
353   }
354   const AliTrackPointArray *GetTrackPointArray() const {
355     return fFriendTrack->GetTrackPointArray(); 
356   }
357   Bool_t RelateToVertexTPC(const AliESDVertex *vtx, Double_t b, Double_t maxd,
358                            AliExternalTrackParam *cParam=0);
359   Bool_t 
360   RelateToVertexTPCBxByBz(const AliESDVertex *vtx, Double_t b[3],Double_t maxd,
361                            AliExternalTrackParam *cParam=0);
362   void GetImpactParametersTPC(Float_t &xy,Float_t &z) const {xy=fdTPC; z=fzTPC;}
363   void GetImpactParametersTPC(Float_t p[2], Float_t cov[3]) const {
364     p[0]=fdTPC; p[1]=fzTPC; cov[0]=fCddTPC; cov[1]=fCdzTPC; cov[2]=fCzzTPC;
365   }
366   Double_t GetConstrainedChi2TPC() const {return fCchi2TPC;}
367
368   Bool_t RelateToVertex(const AliESDVertex *vtx, Double_t b, Double_t maxd,
369                         AliExternalTrackParam *cParam=0);
370   Bool_t 
371   RelateToVertexBxByBz(const AliESDVertex *vtx, Double_t b[3], Double_t maxd,
372                         AliExternalTrackParam *cParam=0);
373   void GetImpactParameters(Float_t &xy,Float_t &z) const {xy=fD; z=fZ;}
374   void GetImpactParameters(Float_t p[2], Float_t cov[3]) const {
375     p[0]=fD; p[1]=fZ; cov[0]=fCdd; cov[1]=fCdz; cov[2]=fCzz;
376   }
377   virtual void Print(Option_t * opt) const ;
378   AliESDEvent* GetESDEvent() const {return fESDEvent;}
379   void         SetESDEvent(AliESDEvent* evt) {fESDEvent = evt;}  
380   //
381   // visualization (M. Ivanov)
382   //
383   void FillPolymarker(TPolyMarker3D *pol, Float_t magf, Float_t minR, Float_t maxR, Float_t stepR);
384
385 protected:
386   
387   AliExternalTrackParam *fCp; // Track parameters constrained to the primary vertex
388   AliExternalTrackParam *fIp; // Track parameters estimated at the inner wall of TPC
389   AliExternalTrackParam *fTPCInner; // Track parameters estimated at the inner wall of TPC using the TPC stand-alone 
390   AliExternalTrackParam *fOp; // Track parameters estimated at the point of maximal radial coordinate reached during the tracking
391   AliExternalTrackParam *fHMPIDp; // Track parameters at HMPID
392   AliESDfriendTrack *fFriendTrack; //! All the complementary information
393
394   TBits    fTPCClusterMap; // Map of clusters, one bit per padrow; 1 if has a cluster on given padrow
395   TBits    fTPCSharedMap;  // Map of clusters, one bit per padrow; 1 if has a shared cluster on given padrow
396
397
398
399   ULong_t   fFlags;          // Reconstruction status flags 
400   Int_t     fID;             // Unique ID of the track
401   Int_t     fLabel;          // Track label
402   Int_t     fITSLabel;       // label according ITS
403   Int_t     fITSModule[12];  // modules crossed by the track in the ITS 
404   Int_t     fTPCLabel;       // label according TPC
405   Int_t     fTRDLabel;       // label according TRD
406   Int_t     fTOFLabel[3];    // TOF label 
407   Int_t     fTOFCalChannel;  // Channel Index of the TOF Signal 
408   Int_t     fTOFindex;       // index of the assigned TOF cluster
409   Int_t     fHMPIDqn;         // 1000000*number of photon clusters + QDC
410   Int_t     fHMPIDcluIdx;     // 1000000*chamber id + cluster idx of the assigned MIP cluster
411   Int_t     fCaloIndex;       // index of associated EMCAL/PHOS cluster (AliESDCaloCluster)
412
413
414   Int_t     fKinkIndexes[3]; // array of indexes of posible kink candidates 
415   Int_t     fV0Indexes[3];   // array of indexes of posible kink candidates 
416
417   Double32_t   fR[AliPID::kSPECIES]; //[0.,0.,8] combined "detector response probability"
418   Double32_t   fITSr[AliPID::kSPECIES]; //[0.,0.,8] "detector response probabilities" (for the PID)
419   Double32_t   fTPCr[AliPID::kSPECIES]; //[0.,0.,8] "detector response probabilities" (for the PID)
420   Double32_t   fTRDr[AliPID::kSPECIES]; //[0.,0.,8] "detector response probabilities" (for the PID)  
421   Double32_t   fTOFr[AliPID::kSPECIES]; //[0.,0.,8] "detector response probabilities" (for the PID)
422   Double32_t   fHMPIDr[AliPID::kSPECIES];//[0.,0.,8] "detector response probabilities" (for the PID)
423
424   Double32_t fHMPIDtrkTheta;//[-2*pi,2*pi,16] theta of the track extrapolated to the HMPID, LORS
425   // how much of this is needed?
426   Double32_t fHMPIDtrkPhi;     //[-2*pi,2*pi,16] phi of the track extrapolated to the HMPID, LORS
427   Double32_t fHMPIDsignal;  // HMPID PID signal (Theta ckov, rad)
428
429   Double32_t   fTrackTime[AliPID::kSPECIES]; // TOFs estimated by the tracking
430   Double32_t   fTrackLength;   // Track length
431
432   Double32_t   fdTPC;          // TPC-only impact parameter in XY plane
433   Double32_t   fzTPC;          // TPC-only impact parameter in Z
434   Double32_t   fCddTPC,fCdzTPC,fCzzTPC; // Covariance matrix of the TPC-only impact parameters 
435   Double32_t   fCchi2TPC;      // [0.,0.,8] TPC-only chi2 at the primary vertex
436
437   Double32_t   fD;             // Impact parameter in XY plane
438   Double32_t   fZ;             // Impact parameter in Z
439   Double32_t   fCdd,fCdz,fCzz; // Covariance matrix of the impact parameters 
440   Double32_t   fCchi2;          // [0.,0.,8] chi2 at the primary vertex
441
442   Double32_t   fITSchi2;        // [0.,0.,8] chi2 in the ITS
443   Double32_t   fTPCchi2;        // [0.,0.,8] chi2 in the TPC
444   Double32_t   fTPCchi2Iter1;  // [0.,0.,8] chi2 in the TPC
445   Double32_t   fTRDchi2;        // [0.,0.,8] chi2 in the TRD
446   Double32_t   fTOFchi2;        // [0.,0.,8] chi2 in the TOF
447   Double32_t fHMPIDchi2;        // [0.,0.,8] chi2 in the HMPID
448
449   Double32_t fGlobalChi2;       // [0.,0.,8] chi2 of the global track
450
451   Double32_t  fITSsignal;     // [0.,0.,10] detector's PID signal
452   Double32_t  fITSdEdxSamples[4]; // [0.,0.,10] ITS dE/dx samples
453
454   Double32_t  fTPCsignal;     // [0.,0.,10] detector's PID signal
455   Double32_t  fTPCsignalS;    // [0.,0.,10] RMS of dEdx measurement
456   Double32_t  fTPCPoints[4];  // [0.,0.,10] TPC points -first, max. dens, last and max density
457
458   Double32_t fTRDsignal;      // detector's PID signal
459   Double32_t fTRDQuality;     // trd quality factor for TOF
460   Double32_t fTRDBudget;      // trd material budget
461
462   Double32_t fTOFsignal;      // detector's PID signal
463   Double32_t fTOFsignalToT;   // detector's ToT signal
464   Double32_t fTOFsignalRaw;   // detector's uncorrected time signal
465   Double32_t fTOFsignalDz;    // local z  of track's impact on the TOF pad 
466   Double32_t fTOFsignalDx;    // local x  of track's impact on the TOF pad 
467   Double32_t fTOFInfo[10];    //! TOF informations
468   Short_t    fTOFdeltaBC;     // detector's Delta Bunch Crossing correction
469   Short_t    fTOFl0l1;        // detector's L0L1 latency correction
470
471   Double32_t fCaloDx ;        // [0.,0.,8] distance to calorimeter cluster in calo plain (phi direction)
472   Double32_t fCaloDz ;        // [0.,0.,8] distance to calorimeter cluster in calo plain (z direction)
473
474   Double32_t fHMPIDtrkX;       // x of the track impact, LORS 
475   Double32_t fHMPIDtrkY;       // y of the track impact, LORS 
476   Double32_t fHMPIDmipX;       // x of the MIP in LORS
477   Double32_t fHMPIDmipY;       // y of the MIP in LORS
478
479
480   UShort_t fTPCncls;       // number of clusters assigned in the TPC
481   UShort_t fTPCnclsF;      // number of findable clusters in the TPC
482   UShort_t fTPCsignalN;    // number of points used for dEdx
483   UShort_t fTPCnclsIter1;  // number of clusters assigned in the TPC - iteration 1
484   UShort_t fTPCnclsFIter1; // number of findable clusters in the TPC - iteration 1
485
486   Char_t  fITSncls;        // number of clusters assigned in the ITS
487   UChar_t fITSClusterMap;  // map of clusters, one bit per a layer
488   UChar_t fTRDncls;        // number of clusters assigned in the TRD
489   UChar_t fTRDncls0;       // number of clusters assigned in the TRD before first material cross
490   UChar_t fTRDntracklets;  // number of TRD tracklets used for tracking/PID
491
492   Int_t fTRDnSlices;     // number of slices used for PID in the TRD
493   Double32_t *fTRDslices;  //[fTRDnSlices] 
494
495   Char_t  fTRDTimBin[kTRDnPlanes];   // Time bin of Max cluster from all six planes
496   Char_t  fVertexID; // ID of the primary vertex this track belongs to
497   AliESDEvent*   fESDEvent; //!Pointer back to event to which the track belongs
498   
499  private:
500
501   AliESDtrack & operator=(const AliESDtrack & );
502   ClassDef(AliESDtrack,57)  //ESDtrack 
503 };
504
505
506
507 #endif 
508