StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StFtpcTrack.hh
1 // $Id: StFtpcTrack.hh,v 1.21 2004/08/05 23:46:30 oldi Exp $
2 // $Log: StFtpcTrack.hh,v $
3 // Revision 1.21 2004/08/05 23:46:30 oldi
4 // The number of points and number of possible points is set by giving the
5 // appropriate detectorId.
6 // A new (inline) function was added to StFtpcTrack, which returns the detectorId.
7 //
8 // Revision 1.20 2004/02/12 19:37:10 oldi
9 // *** empty log message ***
10 //
11 // Revision 1.19 2003/09/16 15:27:02 jcs
12 // removed inline as it would leave a few undefined reference
13 //
14 // Revision 1.18 2003/01/20 09:16:33 oldi
15 // Calculation of residuals simplified.
16 //
17 // Revision 1.17 2003/01/16 18:04:33 oldi
18 // Bugs eliminated. Now it compiles on Solaris again.
19 // Split residuals for global and primary fit.
20 //
21 // Revision 1.16 2002/11/28 09:39:30 oldi
22 // Problem in momentum fit eliminated. Negative vertex Id is not used anymore.
23 // It was used do decide for global or primary fit.
24 // Code was prepared to fill momentum values at outermost points on tracks.
25 // This feature is not used up to now.
26 // Code cleanups.
27 //
28 // Revision 1.15 2002/11/06 13:46:04 oldi
29 // Global/primary fit handling simplified.
30 // Code clean ups.
31 //
32 // Revision 1.14 2002/10/31 13:40:27 oldi
33 // Method GetSector() added.
34 // Method GetMeanR() and GetMeanAlpha() added.
35 //
36 // Revision 1.13 2002/10/24 16:37:45 oldi
37 // dca (impact parameter) is calculated using StHelix::distance(vertexPos), now.
38 // Therefore it is the smallest three dimensional distance of the helix to the
39 // primary vertex.
40 // dca for primary tracks is filled correctly, now. (Even though this value
41 // shouldn't be of great use.)
42 // Code clean-ups.
43 //
44 // Revision 1.12 2002/10/11 15:45:19 oldi
45 // Get FTPC geometry and dimensions from database.
46 // No field fit activated: Returns momentum = 0 but fits a helix.
47 // Bug in TrackMaker fixed (events with z_vertex > outer_ftpc_radius were cut).
48 // QA histograms corrected (0 was supressed).
49 // Code cleanup (several lines of code changed due to *params -> Instance()).
50 // cout -> gMessMgr.
51 //
52 // Revision 1.11 2002/04/29 15:50:05 oldi
53 // All tracking parameters moved to StFtpcTrackingParameters.cc/hh.
54 // In a future version the actual values should be moved to an .idl file (the
55 // interface should be kept as is in StFtpcTrackingParameters.cc/hh).
56 //
57 // Revision 1.10 2002/04/05 16:50:44 oldi
58 // Cleanup of MomentumFit (StFtpcMomentumFit is now part of StFtpcTrack).
59 // Each Track inherits from StHelix, now.
60 // Therefore it is possible to calculate, now:
61 // - residuals
62 // - vertex estimations obtained by back extrapolations of FTPC tracks
63 // Chi2 was fixed.
64 // Many additional minor (and major) changes.
65 //
66 // Revision 1.9 2002/01/29 11:08:16 oldi
67 // Write() renamed to WriteCluster() resp. WriteTrack() to avoid compiler warnings.
68 // As a result the functions TObject::Write() are available again (directly).
69 //
70 // Revision 1.8 2001/01/25 15:22:02 oldi
71 // Review of the complete code.
72 // Fix of several bugs which caused memory leaks:
73 // - Tracks were not allocated properly.
74 // - Tracks (especially split tracks) were not deleted properly.
75 // - TClonesArray seems to have a problem (it could be that I used it in a
76 // wrong way). I changed all occurences to TObjArray which makes the
77 // program slightly slower but much more save (in terms of memory usage).
78 // Speed up of HandleSplitTracks() which is now 12.5 times faster than before.
79 // Cleanup.
80 //
81 // Revision 1.7 2000/11/10 18:38:02 oldi
82 // Changes due to replacement of StThreeVector by TVector3.
83 // Points can be added to a track on either end now.
84 // New data members for dE/dx information.
85 // Cleanup.
86 //
87 // Revision 1.6 2000/07/18 21:22:16 oldi
88 // Changes due to be able to find laser tracks.
89 // Cleanup: - new functions in StFtpcConfMapper, StFtpcTrack, and StFtpcPoint
90 // to bundle often called functions
91 // - short functions inlined
92 // - formulas of StFormulary made static
93 // - avoid streaming of objects of unknown size
94 // (removes the bunch of CINT warnings during compile time)
95 // - two or three minor bugs cured
96 //
97 // Revision 1.5 2000/07/12 11:58:43 jcs
98 // calculate and save FTPC track parameters for unconstrained fit
99 //
100 // Revision 1.4 2000/07/03 12:42:57 jcs
101 // save (pre)Vertex id and unconstrained fit results
102 //
103 // Revision 1.3 2000/06/07 11:43:30 oldi
104 // New data members added: mRowsWithPoints, mChi2Circle, mChi2Length, mRFirst, mRLast, mAlphaFirst, mAlphaLast.
105 // Added the getters and setters for the new data members.
106 // Added GetEta().
107 //
108 // Revision 1.2 2000/05/11 15:14:51 oldi
109 // Changed class names *Hit.* due to already existing class StFtpcHit.cxx in StEvent
110 //
111 // Revision 1.1 2000/05/10 13:39:23 oldi
112 // Initial version of StFtpcTrackMaker
113 //
114 
116 // //
117 // StFtpcTrack class - representation of one FTPC track for the FTPC trackers. //
118 // //
120 
121 #ifndef STAR_StFtpcTrack
122 #define STAR_StFtpcTrack
123 
124 #include "TObjArray.h"
125 #include "TVector3.h"
126 
127 #include "StPhysicalHelix.hh"
128 #include "StThreeVector.hh"
129 #include "StDetectorId.h"
130 
131 #include "MIntArray.h"
132 #include "StFtpcPoint.hh"
133 #include "StFtpcVertex.hh"
134 #include "StFtpcTrackingParams.hh"
135 
136 #include "phys_constants.h"
137 
138 class StFtpcPoint;
139 class StFtpcVertex;
140 
141 class StFtpcTrack : public StPhysicalHelix, public TObject {
142 
143 private:
144 
145  // data from tracker
146  TObjArray *mPoints; // Array of pointers to clusters of track
147  MIntArray *mPointNumbers; // Array of numbers of clusters
148  Int_t mRowsWithPoints; // Binary pattern to know in which row a point is found
149 
150  Int_t mTrackNumber; // number of track
151  Int_t mGlobTrackId; // global track id
152  Double_t mChi2Circle; // Chi squared of circle fit
153  Double_t mChi2Length; // Chi squared of length fit
154  Double_t mTrackLength; // Length of track helix from first to last point
155  Double_t mRadius; // Radius of the helix (projected to a circle)
156  Double_t mCenterX; // x coordinate of the center of the helix (projected to a circle)
157  Double_t mCenterY; // y coordinate of the center of the helix (projected to a circle)
158  Double_t mAlpha0; // angle between first point of track (may be the vertex) and x axis as seen from center of circle projetion
159  Int_t mPid; // particle id
160  Short_t mNMax; // number of possible hits on this track
161  Double_t mRFirst; // radius of (virtual) trackpoint in the first (inner) pad row
162  Double_t mRLast; // radius of (virtual) trackpoint in the last (outer) pad row
163  Double_t mAlphaFirst; // angle of (virtual) trackpoint in the first (inner) pad row
164  Double_t mAlphaLast; // angle of (virtual) trackpoint in the last (outer) pad row
165  Bool_t mFromMainVertex; // true if tracks origin is the main vertex, otherwise false
166 
167  // data from momentum fit
168  TVector3 mP; // ThreeVector of track momentum
169  TVector3 mV; // ThreeVector of vertex used in fit (= first point on track)
170  TVector3 mL; // ThreeVector of last point on track
171  Int_t mQ; // charge measured in fit
172  Double_t mChiSq[2]; // Chi2 of momentum fit
173  Double_t mTheta; // theta value of momentum fit
174  Double_t mDca; // radial impact parameter to main vertex
175 
176  // dE/dx information
177  Double_t mdEdx; // Mean energy loss per length
178  Int_t mNumdEdxHits; // Number of hits accepted for dE/dx
179 
180 
181 public:
182 
183  StFtpcTrack(); // constructor
184  StFtpcTrack(Int_t tracknumber); // constructor which fills tracknumber (all other members are set to default values)
185  virtual ~StFtpcTrack(); // destructor
186 
187  void SetDefaults(); // performs the default setup for the track
188  void AddPoint(StFtpcPoint *point); // adds a point to the track
189  void AddForwardPoint(StFtpcPoint* point); // adds a point after all shifting all existing points by one slot
190  void Fit(); // momentum fit
191  void Fit(StFtpcVertex *vertex, Double_t max_Dca, Bool_t primary_fit); // momentum fit with vertex
192  void CalculateNMax(); // calculates the max. possible number of points
193  Double_t CalcDca(StFtpcVertex *vertex, Bool_t primaryFit); // calculation of distance of closest approach (dca) to main vertex
194  Double_t CalcAlpha0(); // calculation of the angle of xt with respect to the x axis
195  void CalcAndSetAlpha0() { this->SetAlpha0(this->CalcAlpha0()); } // calculates and sets the angle of xt with respect to the x axis
196  void CalcResiduals(Bool_t primaryFit); // calulates the primary or global fit residuals for each point on track
197  void CalcGlobResiduals(); // calulates the global fit residuals for each point on track
198  void CalcPrimResiduals(); // calulates the primary fit residuals for each point on track
199 
200  // momentum fit
201  void MomentumFit(StFtpcVertex *vertex = 0);
202 
203  StThreeVector<Double_t> helixMomentum() const;
204  StThreeVector<Double_t> momentum() const;
205  StThreeVector<Double_t> localMomentum(Double_t s);
206  Double_t chi2Rad() const;
207  Double_t chi2Lin() const;
208 
209  // getter
210  TObjArray *GetHits() const { return mPoints; }
211  MIntArray *GetHitNumbers() const { return mPointNumbers; }
212  Int_t GetRowsWithPoints() const { return mRowsWithPoints; }
213  Int_t GetTrackNumber() const { return mTrackNumber; }
214  Int_t GetGlobalTrackId() const { return mGlobTrackId; }
215  Double_t GetChi2Circle() const { return mChi2Circle; }
216  Double_t GetChi2Length() const { return mChi2Length; }
217  Double_t GetTrackLength() const { return mTrackLength; }
218  Double_t GetRadius() const { return mRadius; }
219  Double_t GetCenterX() const { return mCenterX; }
220  Double_t GetCenterY() const { return mCenterY; }
221  Double_t GetAlpha0() const { return mAlpha0; }
222  Double_t GetPid() const { return mPid; }
223  Short_t GetNMax() const { return mNMax; }
224  Double_t GetRFirst() const { return mRFirst; }
225  Double_t GetRLast() const { return mRLast; }
226  Double_t GetMeanR() const { return TMath::Abs((mRFirst+mRLast)/2.); }
227  Double_t GetAlphaFirst() const { return mAlphaFirst; }
228  Double_t GetAlphaLast() const { return mAlphaLast; }
229  Double_t GetMeanAlpha();
230  Int_t GetNumberOfPoints() const { return mPoints->GetEntriesFast(); }
231  Bool_t ComesFromMainVertex() const { return mFromMainVertex; }
232  TVector3 GetMomentum() const { return mP; }
233  Double_t GetPx() const { return mP.X(); }
234  Double_t GetPy() const { return mP.Y(); }
235  Double_t GetPz() const { return mP.Z(); }
236  Double_t GetP() const;
237  Double_t GetPt() const;
238  Double_t GetPseudoRapidity() const;
239  Double_t GetEta() const;
240  Double_t GetRapidity() const;
241  Int_t GetHemisphere() const;
242  Int_t GetSector() const;
243  StDetectorId GetDetectorId() const { return (GetHemisphere() == 1.) ? kFtpcWestId : kFtpcEastId; }
244 
245  // Note that the counting for the final track is reversed:
246  // The first point is the one closest to the interaction point.
247  TVector3 GetVertex() const { return mV; }
248  TVector3 GetFirstPointOnTrack()const { return mV; }
249  TVector3 GetLastPointOnTrack() const { return mL; }
250  Int_t GetCharge() const { return mQ; }
251  Double_t const *GetChiSq() const { return mChiSq; }
252  Double_t GetTheta() const { return mTheta; }
253  Double_t GetDca() const { return mDca; }
254 
255  Double_t GetdEdx() const { return mdEdx; }
256  Int_t GetNumdEdxHits() const { return mNumdEdxHits; }
257 
258  // setter
259  void SetTrackNumber(Int_t number);
260  void SetGlobalTrackId(Int_t f) { mGlobTrackId = f; }
261  void SetRowsWithPoints(Int_t f) { mRowsWithPoints = f; }
262 
263  void SetPx(Double_t f) { mP.SetX(f); }
264  void SetPy(Double_t f) { mP.SetY(f); }
265  void SetPz(Double_t f) { mP.SetZ(f); }
266 
267  void SetdEdx(Double_t f) { mdEdx = f; }
268  void SetNumdEdxHits(Int_t f) { mNumdEdxHits = f; }
269 
270  void SetChi2Circle(Double_t f) { mChi2Circle = f; }
271  void SetChi2Length(Double_t f) { mChi2Length = f; }
272  void SetRadius(Double_t f) { mRadius = f; }
273  void SetCenterX(Double_t f) { mCenterX = f; }
274  void SetCenterY(Double_t f) { mCenterY = f; }
275  void SetAlpha0(Double_t f) { mAlpha0 = f; }
276  void SetCharge(Int_t f) { mQ = f; }
277  void SetPid(Int_t f) { mPid = f; }
278  void SetRLast(Double_t f) { mRLast = f; }
279  void SetRFirst(Double_t f) { mRFirst = f; }
280  void SetAlphaLast(Double_t f) { mAlphaLast = f; }
281  void SetAlphaFirst(Double_t f) { mAlphaFirst = f; }
282 
283  void SetDca(Double_t f) { mDca = f; }
284  void SetNMax(Short_t f) { mNMax = f; }
285  void ComesFromMainVertex(Bool_t f) { mFromMainVertex = f; }
286 
287  void SetProperties(Bool_t fUsage, Int_t mTrackNumber);
288  void SetPointDependencies();
289 
290 protected:
291 
292  Int_t CircleFit(Double_t x[],Double_t y[], Double_t xw[], Double_t yw[], Int_t num);
293  void LineFit(Double_t *x, Double_t *y, Double_t *z, Double_t *xw, Double_t *yw, Int_t num);
294 
295  StThreeVector<Double_t> mHelixMomentum;
296  StThreeVector<Double_t> mFullMomentum;
297  StThreeVector<Double_t> mVertex;
298  Double_t mXCenter, mYCenter, mFitRadius, mChi2Rad;
299  Double_t mArcOffset, mArcSlope, mChi2Lin;
300  Double_t mZField;
301  Int_t mIterSteps;
302  Int_t mVertexPointOffset;
303 
304  ClassDef(StFtpcTrack, 1) // Ftpc track class
305 };
306 
307 
308 
309 #endif
double x(double s) const
coordinates of helix at point s
Definition: StHelix.hh:182