]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliESDtrack.h
#75811 ZDC: changes to be ported to the release
[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   Double_t GetTPCPoints(Int_t i) const {return fTPCPoints[i];}
217   void    SetKinkIndexes(Int_t points[3]) {
218      for (Int_t i=0;i<3;i++) fKinkIndexes[i] = points[i];
219   }
220   void    SetV0Indexes(Int_t points[3]) {
221      for (Int_t i=0;i<3;i++) fV0Indexes[i] = points[i];
222   }
223   void    SetTPCsignal(Float_t signal, Float_t sigma, UChar_t npoints){ 
224      fTPCsignal = signal; fTPCsignalS = sigma; fTPCsignalN = npoints;
225   }
226   Double_t GetTPCsignal() const {return fTPCsignal;}
227   Double_t GetTPCsignalSigma() const {return fTPCsignalS;}
228   UShort_t GetTPCsignalN() const {return fTPCsignalN;}
229   Double_t GetTPCchi2() const {return fTPCchi2;}
230   Double_t GetTPCchi2Iter1() const {return fTPCchi2Iter1;}
231   UShort_t   GetTPCclusters(Int_t *idx) const;
232   Double_t GetTPCdensity(Int_t row0, Int_t row1) const;
233   Int_t   GetTPCLabel() const {return fTPCLabel;}
234   Int_t   GetKinkIndex(Int_t i) const { return fKinkIndexes[i];}
235   Int_t   GetV0Index(Int_t i) const { return fV0Indexes[i];}
236   const TBits& GetTPCClusterMap() const {return fTPCClusterMap;}
237   const TBits& GetTPCSharedMap() const {return fTPCSharedMap;}
238   void    SetTPCClusterMap(const TBits amap) {fTPCClusterMap = amap;}
239   void    SetTPCSharedMap(const TBits amap) {fTPCSharedMap = amap;}
240
241   void    SetTRDpid(const Double_t *p);
242   
243 // A.Bercuci
244   void    SetTRDntracklets(UChar_t q){fTRDntracklets = q;}
245   UChar_t GetTRDntracklets() const {return (fTRDntracklets>>3)&7;}
246   UChar_t GetTRDntrackletsPID() const {return fTRDntracklets&7;}
247   // TEMPORARY alias asked by the HFE group to allow 
248   // reading of the v4-16-Release data with TRUNK related software (A.Bercuci@Apr 30th 09) 
249   UChar_t GetTRDpidQuality() const {return GetTRDntrackletsPID();}
250 // end A.Bercuci
251
252   void     SetNumberOfTRDslices(Int_t n);
253   Int_t    GetNumberOfTRDslices() const;
254   void     SetTRDslice(Double_t q, Int_t plane, Int_t slice);
255   void     SetTRDmomentum(Double_t p, Int_t plane, Double_t *sp=0x0);
256   Double_t GetTRDslice(Int_t plane, Int_t slice=-1) const;
257   Double_t GetTRDmomentum(Int_t plane, Double_t *sp=0x0) const;
258         
259   void    SetTRDQuality(Float_t quality){fTRDQuality=quality;}
260   Double_t GetTRDQuality()const {return fTRDQuality;}
261   void    SetTRDBudget(Float_t budget){fTRDBudget=budget;}
262   Double_t GetTRDBudget()const {return fTRDBudget;}
263
264   void    SetTRDTimBin(Int_t timbin, Int_t i) {fTRDTimBin[i]=timbin;}
265   void    GetTRDpid(Double_t *p) const;
266   Double_t GetTRDsignal() const {return fTRDsignal;}
267
268   Char_t   GetTRDTimBin(Int_t i) const {return fTRDTimBin[i];}
269   Double_t GetTRDchi2() const {return fTRDchi2;}
270   UChar_t   GetTRDclusters(Int_t *idx) const;
271   UChar_t   GetTRDncls() const {return fTRDncls;}
272   UChar_t   GetTRDncls0() const {return fTRDncls0;}
273   UChar_t   GetTRDtracklets(Int_t *idx) const;
274   void    SetTRDpid(Int_t iSpecies, Float_t p);
275   Double_t GetTRDpid(Int_t iSpecies) const;
276   Int_t   GetTRDLabel() const {return fTRDLabel;}
277
278   void    SetTRDtrack(AliKalmanTrack * track){
279      fFriendTrack->SetTRDtrack(track);
280   }
281   AliKalmanTrack *GetTRDtrack(){
282      return fFriendTrack->GetTRDtrack();
283   }
284
285   void    SetTOFsignal(Double_t tof) {fTOFsignal=tof;}
286   Double_t GetTOFsignal() const {return fTOFsignal;}
287   void    SetTOFsignalToT(Double_t ToT) {fTOFsignalToT=ToT;}
288   Double_t GetTOFsignalToT() const {return fTOFsignalToT;}
289   void    SetTOFsignalRaw(Double_t tof) {fTOFsignalRaw=tof;}
290   Double_t GetTOFsignalRaw() const {return fTOFsignalRaw;}
291   void    SetTOFsignalDz(Double_t dz) {fTOFsignalDz=dz;}
292   Double_t GetTOFsignalDz() const {return fTOFsignalDz;}
293   void    SetTOFsignalDx(Double_t dx) {fTOFsignalDx=dx;}
294   Double_t GetTOFsignalDx() const {return fTOFsignalDx;}
295   void     SetTOFDeltaBC(Short_t deltaBC) {fTOFdeltaBC=deltaBC;};
296   Short_t  GetTOFDeltaBC() const {return fTOFdeltaBC;}
297   void     SetTOFL0L1(Short_t l0l1) {fTOFl0l1=l0l1;};
298   Short_t  GetTOFL0L1() const {return fTOFl0l1;}
299   Double_t GetTOFchi2() const {return fTOFchi2;}
300   void    SetTOFpid(const Double_t *p);
301   void    SetTOFLabel(const Int_t *p);
302   void    GetTOFpid(Double_t *p) const;
303   void    GetTOFLabel(Int_t *p) const;
304   void    GetTOFInfo(Float_t *info) const;
305   void    SetTOFInfo(Float_t *info);
306   Int_t   GetTOFCalChannel() const {return fTOFCalChannel;}
307   Int_t   GetTOFcluster() const {return fTOFindex;}
308   void    SetTOFcluster(Int_t index) {fTOFindex=index;}
309   void    SetTOFCalChannel(Int_t index) {fTOFCalChannel=index;}
310
311 // HMPID methodes +++++++++++++++++++++++++++++++++ (kir)
312   void    SetHMPIDsignal(Double_t theta) {fHMPIDsignal=theta;}
313   Double_t GetHMPIDsignal() const {return fHMPIDsignal;}
314   void    SetHMPIDpid(const Double_t *p);
315   void    GetHMPIDpid(Double_t *p) const;  
316   void    SetHMPIDchi2(Double_t chi2) {fHMPIDchi2=chi2;}
317   Double_t GetHMPIDchi2() const {return fHMPIDchi2;}
318   void    SetHMPIDcluIdx(Int_t ch,Int_t idx) {fHMPIDcluIdx=ch*1000000+idx;}
319   Int_t   GetHMPIDcluIdx() const {return fHMPIDcluIdx;}
320   void    SetHMPIDtrk(Float_t  x, Float_t  y, Float_t  th, Float_t  ph) {
321      fHMPIDtrkX=x; fHMPIDtrkY=y; fHMPIDtrkTheta=th; fHMPIDtrkPhi=ph;
322   }
323   void    GetHMPIDtrk(Float_t &x, Float_t &y, Float_t &th, Float_t &ph) const {
324      x=fHMPIDtrkX; y=fHMPIDtrkY; th=fHMPIDtrkTheta; ph=fHMPIDtrkPhi;
325   }
326   void    SetHMPIDmip(Float_t  x, Float_t  y, Int_t q, Int_t nph=0) {
327      fHMPIDmipX=x; fHMPIDmipY=y; fHMPIDqn=1000000*nph+q;
328   }
329   void    GetHMPIDmip(Float_t &x,Float_t &y,Int_t &q,Int_t &nph) const {
330      x=fHMPIDmipX; y=fHMPIDmipY; q=fHMPIDqn%1000000; nph=fHMPIDqn/1000000;
331   }
332   Bool_t  IsHMPID() const {return fFlags&kHMPIDpid;}
333   Bool_t  IsPureITSStandalone() const {return fFlags&kITSpureSA;}
334   Bool_t  IsMultPrimary() const {return !(fFlags&kMultSec);}
335   Bool_t  IsMultSecondary() const {return (fFlags&kMultSec);}
336
337   Int_t GetEMCALcluster() {return fCaloIndex;}
338   void SetEMCALcluster(Int_t index) {fCaloIndex=index;}
339   Bool_t IsEMCAL() const {return fFlags&kEMCALmatch;}
340
341   Int_t GetPHOScluster() {return fCaloIndex;}
342   void SetPHOScluster(Int_t index) {fCaloIndex=index;}
343   Bool_t IsPHOS() const {return fFlags&kPHOSmatch;}
344   Double_t GetPHOSdx()const{return fCaloDx ;}
345   Double_t GetPHOSdz()const{return fCaloDz ;}
346   void SetPHOSdxdz(Double_t dx, Double_t dz){fCaloDx=dx,fCaloDz=dz;}
347
348
349   void SetTrackPointArray(AliTrackPointArray *points) {
350     fFriendTrack->SetTrackPointArray(points);
351   }
352   const AliTrackPointArray *GetTrackPointArray() const {
353     return fFriendTrack->GetTrackPointArray(); 
354   }
355   Bool_t RelateToVertexTPC(const AliESDVertex *vtx, Double_t b, Double_t maxd,
356                            AliExternalTrackParam *cParam=0);
357   Bool_t 
358   RelateToVertexTPCBxByBz(const AliESDVertex *vtx, Double_t b[3],Double_t maxd,
359                            AliExternalTrackParam *cParam=0);
360   void GetImpactParametersTPC(Float_t &xy,Float_t &z) const {xy=fdTPC; z=fzTPC;}
361   void GetImpactParametersTPC(Float_t p[2], Float_t cov[3]) const {
362     p[0]=fdTPC; p[1]=fzTPC; cov[0]=fCddTPC; cov[1]=fCdzTPC; cov[2]=fCzzTPC;
363   }
364   Double_t GetConstrainedChi2TPC() const {return fCchi2TPC;}
365
366   Bool_t RelateToVertex(const AliESDVertex *vtx, Double_t b, Double_t maxd,
367                         AliExternalTrackParam *cParam=0);
368   Bool_t 
369   RelateToVertexBxByBz(const AliESDVertex *vtx, Double_t b[3], Double_t maxd,
370                         AliExternalTrackParam *cParam=0);
371   void GetImpactParameters(Float_t &xy,Float_t &z) const {xy=fD; z=fZ;}
372   void GetImpactParameters(Float_t p[2], Float_t cov[3]) const {
373     p[0]=fD; p[1]=fZ; cov[0]=fCdd; cov[1]=fCdz; cov[2]=fCzz;
374   }
375   virtual void Print(Option_t * opt) const ;
376   AliESDEvent* GetESDEvent() const {return fESDEvent;}
377   void         SetESDEvent(AliESDEvent* evt) {fESDEvent = evt;}  
378   //
379   // visualization (M. Ivanov)
380   //
381   void FillPolymarker(TPolyMarker3D *pol, Float_t magf, Float_t minR, Float_t maxR, Float_t stepR);
382
383 protected:
384   
385   AliExternalTrackParam *fCp; // Track parameters constrained to the primary vertex
386   AliExternalTrackParam *fIp; // Track parameters estimated at the inner wall of TPC
387   AliExternalTrackParam *fTPCInner; // Track parameters estimated at the inner wall of TPC using the TPC stand-alone 
388   AliExternalTrackParam *fOp; // Track parameters estimated at the point of maximal radial coordinate reached during the tracking
389   AliExternalTrackParam *fHMPIDp; // Track parameters at HMPID
390   AliESDfriendTrack *fFriendTrack; //! All the complementary information
391
392   TBits    fTPCClusterMap; // Map of clusters, one bit per padrow; 1 if has a cluster on given padrow
393   TBits    fTPCSharedMap;  // Map of clusters, one bit per padrow; 1 if has a shared cluster on given padrow
394
395
396
397   ULong_t   fFlags;          // Reconstruction status flags 
398   Int_t     fID;             // Unique ID of the track
399   Int_t     fLabel;          // Track label
400   Int_t     fITSLabel;       // label according ITS
401   Int_t     fITSModule[12];  // modules crossed by the track in the ITS 
402   Int_t     fTPCLabel;       // label according TPC
403   Int_t     fTRDLabel;       // label according TRD
404   Int_t     fTOFLabel[3];    // TOF label 
405   Int_t     fTOFCalChannel;  // Channel Index of the TOF Signal 
406   Int_t     fTOFindex;       // index of the assigned TOF cluster
407   Int_t     fHMPIDqn;         // 1000000*number of photon clusters + QDC
408   Int_t     fHMPIDcluIdx;     // 1000000*chamber id + cluster idx of the assigned MIP cluster
409   Int_t     fCaloIndex;       // index of associated EMCAL/PHOS cluster (AliESDCaloCluster)
410
411
412   Int_t     fKinkIndexes[3]; // array of indexes of posible kink candidates 
413   Int_t     fV0Indexes[3];   // array of indexes of posible kink candidates 
414
415   Double32_t   fR[AliPID::kSPECIES]; //[0.,0.,8] combined "detector response probability"
416   Double32_t   fITSr[AliPID::kSPECIES]; //[0.,0.,8] "detector response probabilities" (for the PID)
417   Double32_t   fTPCr[AliPID::kSPECIES]; //[0.,0.,8] "detector response probabilities" (for the PID)
418   Double32_t   fTRDr[AliPID::kSPECIES]; //[0.,0.,8] "detector response probabilities" (for the PID)  
419   Double32_t   fTOFr[AliPID::kSPECIES]; //[0.,0.,8] "detector response probabilities" (for the PID)
420   Double32_t   fHMPIDr[AliPID::kSPECIES];//[0.,0.,8] "detector response probabilities" (for the PID)
421
422   Double32_t fHMPIDtrkTheta;//[-2*pi,2*pi,16] theta of the track extrapolated to the HMPID, LORS
423   // how much of this is needed?
424   Double32_t fHMPIDtrkPhi;     //[-2*pi,2*pi,16] phi of the track extrapolated to the HMPID, LORS
425   Double32_t fHMPIDsignal;  // HMPID PID signal (Theta ckov, rad)
426
427   Double32_t   fTrackTime[AliPID::kSPECIES]; // TOFs estimated by the tracking
428   Double32_t   fTrackLength;   // Track length
429
430   Double32_t   fdTPC;          // TPC-only impact parameter in XY plane
431   Double32_t   fzTPC;          // TPC-only impact parameter in Z
432   Double32_t   fCddTPC,fCdzTPC,fCzzTPC; // Covariance matrix of the TPC-only impact parameters 
433   Double32_t   fCchi2TPC;      // [0.,0.,8] TPC-only chi2 at the primary vertex
434
435   Double32_t   fD;             // Impact parameter in XY plane
436   Double32_t   fZ;             // Impact parameter in Z
437   Double32_t   fCdd,fCdz,fCzz; // Covariance matrix of the impact parameters 
438   Double32_t   fCchi2;          // [0.,0.,8] chi2 at the primary vertex
439
440   Double32_t   fITSchi2;        // [0.,0.,8] chi2 in the ITS
441   Double32_t   fTPCchi2;        // [0.,0.,8] chi2 in the TPC
442   Double32_t   fTPCchi2Iter1;  // [0.,0.,8] chi2 in the TPC
443   Double32_t   fTRDchi2;        // [0.,0.,8] chi2 in the TRD
444   Double32_t   fTOFchi2;        // [0.,0.,8] chi2 in the TOF
445   Double32_t fHMPIDchi2;        // [0.,0.,8] chi2 in the HMPID
446
447   Double32_t fGlobalChi2;       // [0.,0.,8] chi2 of the global track
448
449   Double32_t  fITSsignal;     // [0.,0.,10] detector's PID signal
450   Double32_t  fITSdEdxSamples[4]; // [0.,0.,10] ITS dE/dx samples
451
452   Double32_t  fTPCsignal;     // [0.,0.,10] detector's PID signal
453   Double32_t  fTPCsignalS;    // [0.,0.,10] RMS of dEdx measurement
454   Double32_t  fTPCPoints[4];  // [0.,0.,10] TPC points -first, max. dens, last and max density
455
456   Double32_t fTRDsignal;      // detector's PID signal
457   Double32_t fTRDQuality;     // trd quality factor for TOF
458   Double32_t fTRDBudget;      // trd material budget
459
460   Double32_t fTOFsignal;      // detector's PID signal
461   Double32_t fTOFsignalToT;   // detector's ToT signal
462   Double32_t fTOFsignalRaw;   // detector's uncorrected time signal
463   Double32_t fTOFsignalDz;    // local z  of track's impact on the TOF pad 
464   Double32_t fTOFsignalDx;    // local x  of track's impact on the TOF pad 
465   Double32_t fTOFInfo[10];    //! TOF informations
466   Short_t    fTOFdeltaBC;     // detector's Delta Bunch Crossing correction
467   Short_t    fTOFl0l1;        // detector's L0L1 latency correction
468
469   Double32_t fCaloDx ;        // [0.,0.,8] distance to calorimeter cluster in calo plain (phi direction)
470   Double32_t fCaloDz ;        // [0.,0.,8] distance to calorimeter cluster in calo plain (z direction)
471
472   Double32_t fHMPIDtrkX;       // x of the track impact, LORS 
473   Double32_t fHMPIDtrkY;       // y of the track impact, LORS 
474   Double32_t fHMPIDmipX;       // x of the MIP in LORS
475   Double32_t fHMPIDmipY;       // y of the MIP in LORS
476
477
478   UShort_t fTPCncls;       // number of clusters assigned in the TPC
479   UShort_t fTPCnclsF;      // number of findable clusters in the TPC
480   UShort_t fTPCsignalN;    // number of points used for dEdx
481   UShort_t fTPCnclsIter1;  // number of clusters assigned in the TPC - iteration 1
482   UShort_t fTPCnclsFIter1; // number of findable clusters in the TPC - iteration 1
483
484   Char_t  fITSncls;        // number of clusters assigned in the ITS
485   UChar_t fITSClusterMap;  // map of clusters, one bit per a layer
486   UChar_t fTRDncls;        // number of clusters assigned in the TRD
487   UChar_t fTRDncls0;       // number of clusters assigned in the TRD before first material cross
488   UChar_t fTRDntracklets;  // number of TRD tracklets used for tracking/PID
489
490   Int_t fTRDnSlices;     // number of slices used for PID in the TRD
491   Double32_t *fTRDslices;  //[fTRDnSlices] 
492
493   Char_t  fTRDTimBin[kTRDnPlanes];   // Time bin of Max cluster from all six planes
494   Char_t  fVertexID; // ID of the primary vertex this track belongs to
495   AliESDEvent*   fESDEvent; //!Pointer back to event to which the track belongs
496   
497  private:
498
499   AliESDtrack & operator=(const AliESDtrack & );
500   ClassDef(AliESDtrack,57)  //ESDtrack 
501 };
502
503
504
505 #endif 
506