Line data Source code
1 : //-*- Mode: C++ -*-
2 : // $Id$
3 : #ifndef ALIHLTGLOBALBARRELTRACK_H
4 : #define ALIHLTGLOBALBARRELTRACK_H
5 : //* This file is property of and copyright by the *
6 : //* ALICE Experiment at CERN, All rights reserved. *
7 : //* See cxx source for full Copyright notice *
8 :
9 : /// @file AliHLTGlobalBarrelTrack.h
10 : /// @author Matthias Richter
11 : /// @date 2009-06-24
12 : /// @brief An AliKalmanTrack implementation for global HLT barrel tracks.
13 : ///
14 :
15 : #include "AliKalmanTrack.h"
16 : #include "AliHLTDataTypes.h"
17 : #include "AliHLTExternalTrackParam.h"
18 : #include <vector>
19 : using std::vector;
20 :
21 : class TClonesArray;
22 : class AliHLTSpacePointContainer;
23 : class AliHLTTrackGeometry;
24 :
25 : /**
26 : * @class AliHLTGlobalBarrelTrack
27 : * Representation of global HLT barrel tracks.
28 : *
29 : * @ingroup alihlt_global_components
30 : */
31 : class AliHLTGlobalBarrelTrack : public AliKalmanTrack
32 : {
33 : public:
34 : /** standard constructor */
35 : AliHLTGlobalBarrelTrack();
36 : /** copy constructor */
37 : AliHLTGlobalBarrelTrack(const AliHLTGlobalBarrelTrack& t);
38 : /** copy constructor */
39 : AliHLTGlobalBarrelTrack(const AliHLTExternalTrackParam& p);
40 : /** copy constructor */
41 : AliHLTGlobalBarrelTrack(const AliExternalTrackParam& p);
42 :
43 : /// assignment operator
44 : /// the standard assignment operator for AliHLTGlobalBarrelTrack is in principle
45 : /// covered by the template definition, however, compiler does not seem to recognize
46 : /// correctly -> effC++ warning,
47 : AliHLTGlobalBarrelTrack& operator=(const AliHLTGlobalBarrelTrack& t) {
48 0 : if (this==&t) return *this;
49 0 : this->~AliHLTGlobalBarrelTrack(); new (this) AliHLTGlobalBarrelTrack(t);
50 0 : return *this;
51 0 : }
52 : template <class c>
53 : AliHLTGlobalBarrelTrack& operator=(const c& t) {
54 0 : this->~AliHLTGlobalBarrelTrack(); new (this) AliHLTGlobalBarrelTrack(t);
55 0 : return *this;
56 0 : }
57 :
58 : /** destructor */
59 : ~AliHLTGlobalBarrelTrack();
60 :
61 : /// inherited from AliKalmanTrack
62 : Int_t GetClusterIndex(Int_t i) const {
63 0 : return (i<(int)fPoints.size()) ?fPoints[i] :0;
64 : }
65 :
66 : /// inherited from AliKalmanTrack, dummy implementation
67 0 : virtual Int_t GetNumberOfTracklets() const {return 0;}
68 : /// inherited from AliKalmanTrack, dummy implementation
69 0 : virtual Int_t GetTrackletIndex(Int_t) const {return -1;}
70 : /// inherited from AliKalmanTrack, dummy implementation
71 0 : virtual Double_t GetPIDsignal() const {return 0.;}
72 :
73 : /// Get the x position of the last assigned point
74 0 : Double_t GetLastPointX() const {return fLastX;}
75 : /// Get the y position of the last assigned point
76 0 : Double_t GetLastPointY() const {return fLastY;}
77 : /// return Track ID
78 0 : Int_t TrackID() const {return fTrackID;}
79 : /// return Track ID, inherited for AliExternalTrackParam
80 0 : Int_t GetID() const {return fTrackID;}
81 :
82 : /// Get the number of associated points
83 : UInt_t GetNumberOfPoints() const;
84 :
85 : /// Get the list of associated points
86 : const UInt_t* GetPoints() const;
87 :
88 : /// Set the space point data
89 0 : void SetSpacePointContainer(AliHLTSpacePointContainer* points) {fSpacePoints=points;}
90 : /// Set track point container
91 : void SetTrackGeometry(AliHLTTrackGeometry* points);
92 : /// Get track point container
93 0 : AliHLTTrackGeometry* GetTrackGeometry() const {return fTrackPoints;}
94 :
95 : /// associate the track space points to the calculated track points
96 : int AssociateSpacePoints(AliHLTTrackGeometry* trackpoints, AliHLTSpacePointContainer& spacepoints) const;
97 :
98 : /// Set the list of associated points
99 : int SetPoints(const UInt_t* pArray, UInt_t arraySize);
100 :
101 : static int ConvertTrackDataArray(const AliHLTTracksData* pTracks, unsigned sizeInByte, vector<AliHLTGlobalBarrelTrack> &tgtArray);
102 :
103 : /// inherited from AliKalmanTrack, dummy implementation
104 0 : Double_t GetPredictedChi2(const AliCluster*) const {return 0.0;}
105 :
106 : /// inherited from AliKalmanTrack, dummy implementation
107 0 : Bool_t PropagateTo(Double_t, Double_t, Double_t) {return kFALSE;}
108 :
109 : /// inherited from AliKalmanTrack, dummy implementation
110 0 : Bool_t Update(const AliCluster*, Double_t, Int_t) {return kFALSE;}
111 :
112 : /// Inherited from TObject, prints the track parameters
113 : virtual void Print(Option_t* option = "") const;
114 :
115 : /// Inherited from TObject, draw the track
116 : virtual void Draw(Option_t *option="");
117 :
118 : int DrawProjXYSpacePoints(Option_t *option, const AliHLTSpacePointContainer* fSpacePoints, const float scale, float center[2]);
119 : int DrawProjXYTrack(Option_t *option, const float scale, float center[2]);
120 : int DrawProjXYTrackPoints(Option_t *option, const float scale, const float center[2], int firstpadrow, int step, float lastpoint[2]);
121 :
122 : Double_t GetPathLengthTo( Double_t x, Double_t b ) const;
123 :
124 : static int ReadTracks(const char* filename, TClonesArray& tgt, AliHLTComponentDataType dt=kAliHLTVoidDataType, unsigned specification=kAliHLTVoidDataSpec);
125 : static int ReadTrackList(const char* listfile, TClonesArray& tgt, AliHLTComponentDataType dt=kAliHLTVoidDataType, unsigned specification=kAliHLTVoidDataSpec);
126 :
127 : /// calculate crossing point with a plane parallel to z axis, at distance x
128 : /// and phi
129 : int CalculateCrossingPoint(float xPlane, float phiPlane, float& u, float& v);
130 :
131 : /// calculate and set internal helix parameters using the global magnetic field
132 : int CalculateHelixParams();
133 : /// calculate and set internal helix parameters
134 : int CalculateHelixParams(float bfield);
135 :
136 : protected:
137 :
138 : private:
139 :
140 : /// array of points
141 : vector<UInt_t> fPoints; //
142 :
143 : /// x position of the last assigned point
144 : Double_t fLastX; //
145 : /// y position of the last assigned point
146 : Double_t fLastY; //
147 :
148 : /// track Id for identification during the reconstruction
149 : Int_t fTrackID; //
150 :
151 : /// helix parameters
152 : Float_t fHelixRadius; //
153 : Float_t fHelixCenterX; //
154 : Float_t fHelixCenterY; //
155 :
156 : /// the space points assigned to the track
157 : AliHLTSpacePointContainer* fSpacePoints; //!
158 : /// the track points according to a geometry
159 : AliHLTTrackGeometry* fTrackPoints; //!
160 :
161 8 : ClassDef(AliHLTGlobalBarrelTrack, 0)
162 : };
163 : #endif
|