]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TPCLib/AliHLTTPCTrack.h
Getting rid of trivial warnings.
[u/mrichter/AliRoot.git] / HLT / TPCLib / AliHLTTPCTrack.h
1 // XEmacs -*-C++-*-
2 // @(#) $Id$
3 // Original: AliHLTTrack.h,v 1.18 2005/03/31 04:48:58 cvetan 
4
5 #ifndef ALIHLTTPCTRACK_H
6 #define ALIHLTTPCTRACK_H
7
8 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
9  * See cxx source for full Copyright notice                               */
10
11 /** @file   AliHLTTPCTrack.h
12     @author Anders Vestbo, Uli Frankenfeld, maintained by Matthias Richter
13     @date   
14     @brief  HLT TPC track base class (conformal mapping)
15 */
16
17 #include "AliTPCtrack.h"
18
19 class AliHLTTPCVertex;
20 class AliHLTTPCSpacePointData;
21
22 /**
23  * @class AliHLTTPCTrack
24  * This class implements the representation of a TPC track, used by the
25  * HLT conformal mapping track finder. <br>
26  * It was originally separated from the offline TPC track class, but in
27  * order to adjust the output format to the offline ESD, AliHLTTPCTrack
28  * now inherits from AliHLTtrack.
29  */
30 class AliHLTTPCTrack : public AliTPCtrack {
31
32  public:
33   
34   AliHLTTPCTrack();
35   virtual ~AliHLTTPCTrack();
36   
37   /**
38    * Copy track parameters.
39    * @param track   pointer to source track
40    */
41   virtual void Copy(AliHLTTPCTrack* track);
42   using AliTPCtrack::Copy;  //TODO: Check if "virtual void Copy(TObject*)" does what it is supposed to do.
43
44   /**
45    * Compare two tracks by the number of hits
46    * @return 0 if equal number of hits, 
47    *         1 if this > track
48    *        -1 if this < track
49    */
50   virtual Int_t Compare(const AliHLTTPCTrack *track) const;
51   using AliTPCtrack::Compare;  //TODO: Check if "virtual Int_t Compare(TObject*)" does what it is supposed to do.
52
53   /**
54    * Fit the assigned spacepoints to a helix.
55    * The function sets teh track parameters.
56    */
57   virtual void CalculateHelix();
58   
59   Bool_t CalculateReferencePoint(Double_t angle,Double_t radius=132);//Calculate Reference Point
60   Bool_t CalculateEdgePoint(Double_t angle);//Calculate crossing point with line
61   Bool_t CalculatePoint(Double_t xplane);   //Calculate crossing point with X-plane
62   Bool_t IsPoint() {return fIsPoint;}
63   Double_t GetCrossingAngle(Int_t padrow,Int_t slice=-1);
64   Bool_t GetCrossingPoint(Int_t padrow,Float_t *xyz);
65   Double_t GetDistance(Double_t /*x0*/,Double_t /*x1*/){return 0;}
66   void UpdateToFirstPoint();
67   void GetClosestPoint(AliHLTTPCVertex *vertex,Double_t &closestX,Double_t &closestY,Double_t &closestZ);
68   void Rotate(Int_t slice,Bool_t tolocal=kFALSE);
69   Bool_t IsLocal() const {return fIsLocal;}
70   virtual void Print(Option_t* option = "") const;
71   using AliTPCtrack::Print;
72
73   // getter
74   Double_t GetFirstPointX() const {return fFirstPoint[0];}
75   Double_t GetFirstPointY() const {return fFirstPoint[1];}
76   Double_t GetFirstPointZ() const {return fFirstPoint[2];}
77   Double_t GetLastPointX() const {return fLastPoint[0];}
78   Double_t GetLastPointY() const {return fLastPoint[1];}
79   Double_t GetLastPointZ() const {return fLastPoint[2];}
80
81   Double_t GetPointPsi() const {return fPointPsi;}
82   Double_t GetPointX() const {return fPoint[0];}
83   Double_t GetPointY() const {return fPoint[1];}
84   Double_t GetPointZ() const {return fPoint[2];}
85
86   Double_t GetPt() const {return fPt;}
87   Double_t GetTgl() const {return fTanl;}
88   Double_t GetPsi() const {return fPsi;}
89   Double_t GetPhi0() const {return fPhi0;}
90   Double_t GetR0() const {return fR0;}
91   Double_t GetZ0() const {return fFirstPoint[2];}
92   Float_t GetPID() const {return fPID;}
93
94   Double_t GetPterr() const {return fPterr;}
95   Double_t GetPsierr() const {return fPsierr;}
96   Double_t GetTglerr() const {return fTanlerr;}
97   
98   Double_t GetKappa() const {return fKappa;}
99   Double_t GetRadius() const {return fRadius;}
100   Double_t GetCenterX() const {return fCenterX;}
101   Double_t GetCenterY() const {return fCenterY;}
102
103   Int_t GetNHits() const {return fNHits;}
104   Int_t   GetNumberOfPoints()   const {return fNHits;}
105   Bool_t  ComesFromMainVertex() const {return fFromMainVertex;}
106     
107   Double_t GetPx() const {return fPt*cos(fPsi);}
108   Double_t GetPy() const {return fPt*sin(fPsi);}
109   Double_t GetPz() const {return fPt*fTanl;}
110   
111   Double_t GetP() const;
112   Double_t GetPseudoRapidity() const;
113   Double_t GetRapidity() const;
114   
115   Int_t GetCharge() const {return fQ;}
116   Int_t GetMCid() const {return fMCid;}
117   Double_t GetLength() const {return fLength;}
118
119   Int_t GetFirstRow() const {return fRowRange[0];}
120   Int_t GetLastRow()  const {return fRowRange[1];}
121   Int_t GetSector()   const {return fSector;}
122
123   UInt_t *GetHitNumbers() {return fHitNumbers;}
124
125   // setter   
126   void SetPID(Float_t pid) {fPID=pid;}  
127   void SetMCid(Int_t f) {fMCid = f;}
128   void SetFirstPoint(Double_t f,Double_t g,Double_t h) {fFirstPoint[0]=f; fFirstPoint[1]=g; fFirstPoint[2]=h;}
129   void SetLastPoint(Double_t f,Double_t g,Double_t h) {fLastPoint[0]=f; fLastPoint[1]=g; fLastPoint[2]=h;}
130   void SetHits(Int_t nhits,UInt_t *hits) {memcpy(fHitNumbers,hits,nhits*sizeof(UInt_t));}
131   void SetPhi0(Double_t f) {fPhi0 = f;}
132   void SetPsi(Double_t f) {fPsi = f;}
133   void SetR0(Double_t f) {fR0 = f;}
134   void SetTgl(Double_t f) {fTanl =f;}
135   void SetZ0(Double_t f) {fFirstPoint[2] = f;}
136   void SetPt(Double_t f) {fPt = f;}
137   void SetLength(Double_t f) {fLength = f;}
138   void SetPterr(Double_t f) {fPterr = f;}
139   void SetPsierr(Double_t f) {fPsierr = f;}
140   void SetZ0err(Double_t f) {fZ0err = f;}
141   void SetTglerr(Double_t f) {fTanlerr = f;}
142   void SetKappa(Double_t f) {fKappa = f;}
143   void SetNHits(Int_t f) {fNHits = f;}
144   void SetRowRange(Int_t f,Int_t g) {fRowRange[0]=f; fRowRange[1]=g;}
145   void SetSector(Int_t f) {fSector = f;}
146   void SetRadius(Double_t f) {fRadius = f;}
147   void SetCenterX(Double_t f) {fCenterX = f;}
148   void SetCenterY(Double_t f) {fCenterY = f;}
149   void SetCharge(Int_t f) {fQ = f;}
150   void ComesFromMainVertex(Bool_t f) {fFromMainVertex = f;}
151
152   /**
153    * Convert all track parameters to the format of AliKalmanTrack
154    * The AliKalmanTrack class implements the track parametrization for offline ITS, TPC
155    * and TRD tracking. The function calculates and sets the parameters of the
156    * parent class (Note: AliHLTTPCTrack inherits from AliTPCtrack and thus
157    * AliKalmanTrack).
158    */
159   int Convert2AliKalmanTrack();
160
161  private:
162
163   Int_t fNHits; //Number of hits
164   Int_t fMCid;  //Assigned id from MC data.
165
166   Double_t fKappa;   // Signed curvature (projected to a circle)
167   Double_t fRadius;  // Radius of the helix (projected to a circle)
168   Double_t fCenterX; // x coordinate of the center of the helix (projected to a circle)
169   Double_t fCenterY; // y coordinate of the center of the helix (projected to a circle)
170   Bool_t   fFromMainVertex; // true if tracks origin is the main vertex, otherwise false
171   
172   Int_t fRowRange[2]; //Subsector where this track was build
173   Int_t fSector;      //Sector # where  this track was build
174
175   //data from momentum fit
176   Int_t    fQ;    //charge measured fit
177     
178   //track parameters:
179   Double_t fTanl; //tan of dipangle
180   Double_t fPsi;  //azimuthal angle of the momentum 
181   Double_t fPt;   //transverse momentum
182   Double_t fLength; //length of track (s)
183   
184   Double_t fPterr;   //error in pt
185   Double_t fPsierr;  //error in psi
186   Double_t fZ0err;   //error in first point
187   Double_t fTanlerr; //error in tanl
188
189   Double_t fPhi0; //azimuthal angle of the first point
190   Double_t fR0;   //radius of the first point
191   Double_t fZ0;   //z coordinate of the first point (fFirstPoint[2])
192
193   Double_t fFirstPoint[3]; //first point
194   Double_t fLastPoint[3];  //last point
195   Double_t fPoint[3]; //point
196   Double_t fPointPsi; //azimuthal angle of the momentum at Point
197
198   Bool_t fIsPoint;  //Helix crosses the X-plane
199   Bool_t fIsLocal; //Track given in local coordinates.
200
201   Float_t fPID; //pid 
202   UInt_t fHitNumbers[159]; //Array of hit numbers for this track
203
204   Bool_t IsPoint(Bool_t ispoint) {fIsPoint = ispoint;return fIsPoint;}
205   
206   ClassDef(AliHLTTPCTrack,2) //Base track class
207 };
208 #endif