QwAnalysis
QwTrack.h
Go to the documentation of this file.
1 /*!
2  * \file QwTrack.h
3  * \brief Definition of the track class
4  *
5  * \author Wouter Deconinck
6  * \date 2009-12-11
7  */
8 
9 #ifndef QWTRACK_H
10 #define QWTRACK_H
11 
12 // Qweak headers
13 #include "VQwTrackingElement.h"
14 #include "QwPartialTrack.h"
15 #include "QwVertex.h"
16 
17 /**
18  * \class QwTrack
19  * \ingroup QwTracking
20  * \brief Contains the complete track as a concatenation of partial tracks
21  *
22  * A QwTrack contains the complete description of the track as a concatenation
23  * of multiple partial tracks. Associated with the track is the kinematical
24  * information, for use in the final Q^2 determination.
25  *
26  * Several vectors of QwPartialTracks are stored in the QwTrack object. This
27  * allows for combining different QwPartialTracks with each other, and selecting
28  * the optimal fit.
29  */
30 class QwTrack: public VQwTrackingElement, public QwObjectCounter<QwTrack> {
31 
32  private:
33 
34  //! Number of partial tracks in this track
36  //! List of partial tracks in this track
37  std::vector<QwPartialTrack*> fQwPartialTracks;
38 
39  public:
40 
41  /// Default constructor
42  QwTrack();
43  /// Constructor with front and back partial track
44  QwTrack(const QwPartialTrack* front, const QwPartialTrack* back);
45  /// Copy constructor by reference
46  QwTrack(const QwTrack& that);
47  /// Copy constructor from pointer
48  QwTrack(const QwTrack* that);
49  /// Virtual destructor
50  virtual ~QwTrack();
51 
52  /// Assignment operator
53  QwTrack& operator=(const QwTrack& that);
54 
55  /// Initialization
56  void Initialize();
57 
58  ///@{
59  /// Creating and adding partial tracks
61  void AddPartialTrack(const QwPartialTrack* partialtrack);
62  void AddPartialTrackList(const std::vector<QwPartialTrack*> &partialtracklist);
63  void ClearPartialTracks(Option_t *option = "");
64  void ResetPartialTracks(Option_t *option = "");
65  ///@}
66  // Get the number of partial tracks
67  Int_t GetNumberOfPartialTracks() const { return fNQwPartialTracks; };
68  //! \brief Get the specified partial track
69  const QwPartialTrack* GetPartialTrack(const int pt) const { return fQwPartialTracks.at(pt); };
70 
71  // Print the list of partial tracks
72  void PrintPartialTracks(Option_t *option = "") const;
73 
74 
75  /// Set magnetic field integral
76  void SetMagneticFieldIntegral(double bdl, double bdlx, double bdly, double bdlz) {
77  fBdl = bdl;
78  fBdlXYZ = TVector3(bdlx,bdly,bdlz);
79  }
80  /// Get magnetic field integral
81  const TVector3& GetMagneticFieldIntegral() const {
82  return fBdlXYZ;
83  }
84 
85 
86  /// Output stream operator for tracks
87  friend std::ostream& operator<< (std::ostream& stream, const QwTrack& t);
88 
89  public:
90 
91  //@{
92  /// Quantities determined from front partial track
93  double fPhi; ///< Azimuthal angle phi of track at primary vertex
94  double fTheta; ///< Polar angle theta of track at primary vertex
95  double fVertexZ; ///< Primary vertex position in longitudinal direction
96  double fVertexR; ///< Primary vertex position in transverse direction
97  //@}
98 
99  double fChi; ///< Combined chi, i.e. the sum of the chi for front and back partial track
100  double fMomentum; ///< Spectrometer momentum
101 
102  int fIterationsNewton; ///< Number of iterations in Newton's method
103  int fIterationsRungeKutta; ///< Number of iterations in Runge-Kutta method
104 
105  //@{
106  /// Matching of front and back track position and direction at matching plane
107  TVector3 fPositionDiff; ///< Difference in position at matching plane
108  double fPositionXoff; ///< Difference in X position at matching plane
109  double fPositionYoff; ///< Difference in Y position at matching plane
110  double fPositionRoff; ///< Difference in radial position at matching plane
111  double fPositionPhioff; ///< Difference in azimuthal angle phi at matching plane
112  double fPositionThetaoff; ///< Difference in polar angle theta at matching plane
113  TVector3 fDirectionDiff; ///< Difference in momentum at matching plane
114  double fDirectionXoff; ///< Difference in X momentum at matching plane
115  double fDirectionYoff; ///< Difference in Y momentum at matching plane
116  double fDirectionZoff; ///< Difference in Z momentum at matching plane
117  double fDirectionPhioff; ///< Difference in momentum azimuthal angle at matching plane
118  double fDirectionThetaoff; ///< Difference in momentum polar angle at matching plane
119  //@}
120 
121  //@{
122  /// Position and direction before and after swimming
123  TVector3 fStartPosition; ///< Start position of front track
124  TVector3 fStartDirection; ///< Start direction of front track
125 
126  TVector3 fEndPositionGoal; ///< Goal position of back track
127  TVector3 fEndDirectionGoal; ///< Goal direction of back track
128 
129  TVector3 fEndPositionActual; ///< Actual position of track at back plane
130  TVector3 fEndDirectionActual; ///< Actual direction of track at back plane
131 
132  int fIterationsRK4; ///< Number of iterations using Runge-Kutta 4th order
133  TVector3 fEndPositionActualRK4; ///< Actual position of track at back plane using Runge-Kutta 4th order
134  TVector3 fEndDirectionActualRK4; ///< Actual direction of track at back plane using Runge-Kutta 4th order
135 
136  int fIterationsRKF45; ///< Number of iterations using Runge-Kutta-Fehlberg
137  TVector3 fEndPositionActualRKF45; ///< Actual position of track at back plane using Runge-Kutta-Fehlberg
138  TVector3 fEndDirectionActualRKF45; ///< Actual direction of track at back plane using Runge-Kutta-Fehlberg
139  //@}
140 
141  //@{
142  /// Magnetic field integrals
143  double fBdl; ///< Magnetic field integral along track
144  TVector3 fBdlXYZ; ///< Magnetic field integral along track
145  //@}
146 
147  QwPartialTrack *fFront; ///< Front partial track
148  QwPartialTrack *fBack; ///< Back partial track
149 
150  ClassDef(QwTrack,3);
151 
152 }; // class QwTrack
153 
154 
155 /**
156  * Method to print vectors conveniently
157  * @param stream Output stream
158  * @param v Vector
159  * @return Output stream
160  */
161 inline std::ostream& operator<< (std::ostream& stream, const TVector3& v)
162 {
163  stream << "(" << v.X() << "," << v.Y() << "," << v.Z() << ")";
164  return stream;
165 }
166 
167 
168 #endif // QWTRACK_H
double fDirectionYoff
Difference in Y momentum at matching plane.
Definition: QwTrack.h:115
double fChi
Combined chi, i.e. the sum of the chi for front and back partial track.
Definition: QwTrack.h:99
double fBdl
Magnetic field integrals.
Definition: QwTrack.h:143
int fIterationsRungeKutta
Number of iterations in Runge-Kutta method.
Definition: QwTrack.h:103
double fDirectionThetaoff
Difference in momentum polar angle at matching plane.
Definition: QwTrack.h:118
std::ostream & operator<<(std::ostream &out, const QwColor &color)
Output stream operator which uses the enum-to-escape-code mapping.
Definition: QwColor.h:153
double fPositionPhioff
Difference in azimuthal angle phi at matching plane.
Definition: QwTrack.h:111
int fIterationsRK4
Number of iterations using Runge-Kutta 4th order.
Definition: QwTrack.h:132
int fIterationsNewton
Number of iterations in Newton&#39;s method.
Definition: QwTrack.h:102
TVector3 fBdlXYZ
Magnetic field integral along track.
Definition: QwTrack.h:144
QwTrack & operator=(const QwTrack &that)
Assignment operator.
Definition: QwTrack.cc:91
TVector3 fStartPosition
Position and direction before and after swimming.
Definition: QwTrack.h:123
TVector3 fEndPositionActualRK4
Actual position of track at back plane using Runge-Kutta 4th order.
Definition: QwTrack.h:133
TVector3 fDirectionDiff
Difference in momentum at matching plane.
Definition: QwTrack.h:113
Definition of the partial track class.
double fDirectionPhioff
Difference in momentum azimuthal angle at matching plane.
Definition: QwTrack.h:117
Memory management structure to count objects.
double fPositionRoff
Difference in radial position at matching plane.
Definition: QwTrack.h:110
virtual ~QwTrack()
Virtual destructor.
Definition: QwTrack.cc:78
void SetMagneticFieldIntegral(double bdl, double bdlx, double bdly, double bdlz)
Set magnetic field integral.
Definition: QwTrack.h:76
TVector3 fEndDirectionActual
Actual direction of track at back plane.
Definition: QwTrack.h:130
TVector3 fStartDirection
Start direction of front track.
Definition: QwTrack.h:124
const TVector3 & GetMagneticFieldIntegral() const
Get magnetic field integral.
Definition: QwTrack.h:81
void ResetPartialTracks(Option_t *option="")
Definition: QwTrack.cc:197
QwTrack()
Default constructor.
Definition: QwTrack.cc:7
Int_t GetNumberOfPartialTracks() const
Definition: QwTrack.h:67
friend std::ostream & operator<<(std::ostream &stream, const QwTrack &t)
Output stream operator for tracks.
Definition: QwTrack.cc:215
Int_t fNQwPartialTracks
Number of partial tracks in this track.
Definition: QwTrack.h:35
ClassDef(QwTrack, 3)
TVector3 fEndPositionActual
Actual position of track at back plane.
Definition: QwTrack.h:129
Contains the complete track as a concatenation of partial tracks.
Definition: QwTrack.h:30
double fDirectionXoff
Difference in X momentum at matching plane.
Definition: QwTrack.h:114
QwPartialTrack * CreateNewPartialTrack()
Definition: QwTrack.cc:163
double fDirectionZoff
Difference in Z momentum at matching plane.
Definition: QwTrack.h:116
TVector3 fEndPositionActualRKF45
Actual position of track at back plane using Runge-Kutta-Fehlberg.
Definition: QwTrack.h:137
void PrintPartialTracks(Option_t *option="") const
Definition: QwTrack.cc:203
void Initialize()
Initialization.
Definition: QwTrack.cc:149
TVector3 fPositionDiff
Matching of front and back track position and direction at matching plane.
Definition: QwTrack.h:107
void AddPartialTrack(const QwPartialTrack *partialtrack)
Definition: QwTrack.cc:171
void AddPartialTrackList(const std::vector< QwPartialTrack * > &partialtracklist)
Definition: QwTrack.cc:178
TVector3 fEndDirectionGoal
Goal direction of back track.
Definition: QwTrack.h:127
double fPhi
Quantities determined from front partial track.
Definition: QwTrack.h:93
QwPartialTrack * fBack
Back partial track.
Definition: QwTrack.h:148
int fIterationsRKF45
Number of iterations using Runge-Kutta-Fehlberg.
Definition: QwTrack.h:136
TVector3 fEndPositionGoal
Goal position of back track.
Definition: QwTrack.h:126
double fPositionXoff
Difference in X position at matching plane.
Definition: QwTrack.h:108
Definition of virtual base class for all tracking elements.
QwPartialTrack * fFront
Front partial track.
Definition: QwTrack.h:147
TVector3 fEndDirectionActualRK4
Actual direction of track at back plane using Runge-Kutta 4th order.
Definition: QwTrack.h:134
double fPositionYoff
Difference in Y position at matching plane.
Definition: QwTrack.h:109
TVector3 fEndDirectionActualRKF45
Actual direction of track at back plane using Runge-Kutta-Fehlberg.
Definition: QwTrack.h:138
double fPositionThetaoff
Difference in polar angle theta at matching plane.
Definition: QwTrack.h:112
double fTheta
Polar angle theta of track at primary vertex.
Definition: QwTrack.h:94
std::vector< QwPartialTrack * > fQwPartialTracks
List of partial tracks in this track.
Definition: QwTrack.h:37
double fMomentum
Spectrometer momentum.
Definition: QwTrack.h:100
const QwPartialTrack * GetPartialTrack(const int pt) const
Get the specified partial track.
Definition: QwTrack.h:69
double fVertexZ
Primary vertex position in longitudinal direction.
Definition: QwTrack.h:95
Contains the straight part of a track in one region only.
void ClearPartialTracks(Option_t *option="")
Definition: QwTrack.cc:186
double fVertexR
Primary vertex position in transverse direction.
Definition: QwTrack.h:96
Virtual base class for all tracking elements.