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