QwAnalysis
QwTrack Class Reference

Contains the complete track as a concatenation of partial tracks. More...

#include <QwTrack.h>

+ Inheritance diagram for QwTrack:
+ Collaboration diagram for QwTrack:

Public Member Functions

 QwTrack ()
 Default constructor. More...
 
 QwTrack (const QwPartialTrack *front, const QwPartialTrack *back)
 Constructor with front and back partial track. More...
 
 QwTrack (const QwTrack &that)
 Copy constructor by reference. More...
 
 QwTrack (const QwTrack *that)
 Copy constructor from pointer. More...
 
virtual ~QwTrack ()
 Virtual destructor. More...
 
QwTrackoperator= (const QwTrack &that)
 Assignment operator. More...
 
void Initialize ()
 Initialization. More...
 
Int_t GetNumberOfPartialTracks () const
 
const QwPartialTrackGetPartialTrack (const int pt) const
 Get the specified partial track. More...
 
void PrintPartialTracks (Option_t *option="") const
 
void SetMagneticFieldIntegral (double bdl, double bdlx, double bdly, double bdlz)
 Set magnetic field integral. More...
 
const TVector3 & GetMagneticFieldIntegral () const
 Get magnetic field integral. More...
 
 ClassDef (QwTrack, 3)
 
QwPartialTrackCreateNewPartialTrack ()
 
void AddPartialTrack (const QwPartialTrack *partialtrack)
 
void AddPartialTrackList (const std::vector< QwPartialTrack * > &partialtracklist)
 
void ClearPartialTracks (Option_t *option="")
 
void ResetPartialTracks (Option_t *option="")
 
- Public Member Functions inherited from VQwTrackingElement
 VQwTrackingElement ()
 Default constructor. More...
 
 VQwTrackingElement (const VQwTrackingElement &that)
 
virtual ~VQwTrackingElement ()
 Virtual destructor. More...
 
VQwTrackingElementoperator= (const VQwTrackingElement &that)
 Assignment operator. More...
 
const QwDetectorInfoGetDetectorInfo () const
 Get the detector info pointer. More...
 
void SetDetectorInfo (const QwDetectorInfo *detectorinfo)
 Set the detector info pointer. More...
 
EQwRegionID GetRegion () const
 Get the region. More...
 
void SetRegion (EQwRegionID region)
 Set the region. More...
 
EQwDetectorPackage GetPackage () const
 Get the package. More...
 
void SetPackage (EQwDetectorPackage package)
 Set the package. More...
 
int GetOctant () const
 Get the octant number. More...
 
void SetOctant (int octant)
 Set the octant number. More...
 
EQwDirectionID GetDirection () const
 Get the direction. More...
 
void SetDirection (EQwDirectionID direction)
 Set the direction. More...
 
int GetPlane () const
 Get the plane number. More...
 
void SetPlane (int plane)
 Set the plane number. More...
 
int GetElement () const
 Get the element number. More...
 
void SetElement (int element)
 Set the element number. More...
 
void SetGeometryTo (const VQwTrackingElement &e)
 Copy the geometry info from another object. More...
 
- Public Member Functions inherited from QwObjectCounter< QwTrack >
 QwObjectCounter ()
 Default constructor. More...
 
 QwObjectCounter (const QwObjectCounter &)
 Copy constructor. More...
 
virtual ~QwObjectCounter ()
 Destructor. More...
 

Data Fields

double fChi
 Combined chi, i.e. the sum of the chi for front and back partial track. More...
 
double fMomentum
 Spectrometer momentum. More...
 
int fIterationsNewton
 Number of iterations in Newton's method. More...
 
int fIterationsRungeKutta
 Number of iterations in Runge-Kutta method. More...
 
QwPartialTrackfFront
 Front partial track. More...
 
QwPartialTrackfBack
 Back partial track. More...
 
double fPhi
 Quantities determined from front partial track. More...
 
double fTheta
 Polar angle theta of track at primary vertex. More...
 
double fVertexZ
 Primary vertex position in longitudinal direction. More...
 
double fVertexR
 Primary vertex position in transverse direction. More...
 
TVector3 fPositionDiff
 Matching of front and back track position and direction at matching plane. More...
 
double fPositionXoff
 Difference in X position at matching plane. More...
 
double fPositionYoff
 Difference in Y position at matching plane. More...
 
double fPositionRoff
 Difference in radial position at matching plane. More...
 
double fPositionPhioff
 Difference in azimuthal angle phi at matching plane. More...
 
double fPositionThetaoff
 Difference in polar angle theta at matching plane. More...
 
TVector3 fDirectionDiff
 Difference in momentum at matching plane. More...
 
double fDirectionXoff
 Difference in X momentum at matching plane. More...
 
double fDirectionYoff
 Difference in Y momentum at matching plane. More...
 
double fDirectionZoff
 Difference in Z momentum at matching plane. More...
 
double fDirectionPhioff
 Difference in momentum azimuthal angle at matching plane. More...
 
double fDirectionThetaoff
 Difference in momentum polar angle at matching plane. More...
 
TVector3 fStartPosition
 Position and direction before and after swimming. More...
 
TVector3 fStartDirection
 Start direction of front track. More...
 
TVector3 fEndPositionGoal
 Goal position of back track. More...
 
TVector3 fEndDirectionGoal
 Goal direction of back track. More...
 
TVector3 fEndPositionActual
 Actual position of track at back plane. More...
 
TVector3 fEndDirectionActual
 Actual direction of track at back plane. More...
 
int fIterationsRK4
 Number of iterations using Runge-Kutta 4th order. More...
 
TVector3 fEndPositionActualRK4
 Actual position of track at back plane using Runge-Kutta 4th order. More...
 
TVector3 fEndDirectionActualRK4
 Actual direction of track at back plane using Runge-Kutta 4th order. More...
 
int fIterationsRKF45
 Number of iterations using Runge-Kutta-Fehlberg. More...
 
TVector3 fEndPositionActualRKF45
 Actual position of track at back plane using Runge-Kutta-Fehlberg. More...
 
TVector3 fEndDirectionActualRKF45
 Actual direction of track at back plane using Runge-Kutta-Fehlberg. More...
 
double fBdl
 Magnetic field integrals. More...
 
TVector3 fBdlXYZ
 Magnetic field integral along track. More...
 

Private Attributes

Int_t fNQwPartialTracks
 Number of partial tracks in this track. More...
 
std::vector< QwPartialTrack * > fQwPartialTracks
 List of partial tracks in this track. More...
 

Friends

std::ostream & operator<< (std::ostream &stream, const QwTrack &t)
 Output stream operator for tracks. More...
 

Additional Inherited Members

- Static Public Member Functions inherited from QwObjectCounter< QwTrack >
static size_t GetObjectsCreated ()
 Get number of objects ever created. More...
 
static size_t GetObjectsAlive ()
 Get number of objects still alive. More...
 
- Protected Member Functions inherited from VQwTrackingElement
 ClassDef (VQwTrackingElement, 1)
 
- Protected Attributes inherited from VQwTrackingElement
const QwDetectorInfofDetectorInfo
 
EQwRegionID fRegion
 ///< Detector info pointer More...
 
EQwDetectorPackage fPackage
 Package. More...
 
int fOctant
 Octant number. More...
 
EQwDirectionID fDirection
 Direction. More...
 
int fPlane
 Plane number. More...
 
int fElement
 Element number. More...
 

Detailed Description

Contains the complete track as a concatenation of partial tracks.

A QwTrack contains the complete description of the track as a concatenation of multiple partial tracks. Associated with the track is the kinematical information, for use in the final Q^2 determination.

Several vectors of QwPartialTracks are stored in the QwTrack object. This allows for combining different QwPartialTracks with each other, and selecting the optimal fit.

Definition at line 30 of file QwTrack.h.

Constructor & Destructor Documentation

QwTrack::QwTrack ( )

Default constructor.

Default constructor

Definition at line 7 of file QwTrack.cc.

References Initialize().

8 {
9  // Initialize
10  Initialize();
11 }
void Initialize()
Initialization.
Definition: QwTrack.cc:149

+ Here is the call graph for this function:

QwTrack::QwTrack ( const QwPartialTrack front,
const QwPartialTrack back 
)

Constructor with front and back partial track.

Constructor with front and back partial track

Parameters
frontFront partial track
backBack partial track

Definition at line 18 of file QwTrack.cc.

References AddPartialTrack(), fBack, fFront, VQwTrackingElement::GetPackage(), Initialize(), and VQwTrackingElement::SetPackage().

19 {
20  // Initialize
21  Initialize();
22 
23  // Null pointer
24  if (front == 0 || back == 0) return;
25 
26  if (front->GetPackage() == back->GetPackage())
27  SetPackage(front->GetPackage());
28 
29  // Copy tracks
30  fFront = new QwPartialTrack(front);
31  fBack = new QwPartialTrack(back);
32 
33  // Add partial tracks
34  AddPartialTrack(front);
35  AddPartialTrack(back);
36 }
void SetPackage(EQwDetectorPackage package)
Set the package.
void Initialize()
Initialization.
Definition: QwTrack.cc:149
void AddPartialTrack(const QwPartialTrack *partialtrack)
Definition: QwTrack.cc:171
QwPartialTrack * fBack
Back partial track.
Definition: QwTrack.h:148
QwPartialTrack * fFront
Front partial track.
Definition: QwTrack.h:147
Contains the straight part of a track in one region only.
EQwDetectorPackage GetPackage() const
Get the package.

+ Here is the call graph for this function:

QwTrack::QwTrack ( const QwTrack that)

Copy constructor by reference.

Copy constructor by reference

Parameters
thatOriginal track

Definition at line 42 of file QwTrack.cc.

References fBack, fFront, and Initialize().

43 : VQwTrackingElement(that)
44 {
45  // Initialize
46  Initialize();
47 
48  *this = that;
49 
50  // Copy tracks
51  fFront = new QwPartialTrack(that.fFront);
52  fBack = new QwPartialTrack(that.fBack);
53 }
void Initialize()
Initialization.
Definition: QwTrack.cc:149
VQwTrackingElement()
Default constructor.
QwPartialTrack * fBack
Back partial track.
Definition: QwTrack.h:148
QwPartialTrack * fFront
Front partial track.
Definition: QwTrack.h:147
Contains the straight part of a track in one region only.

+ Here is the call graph for this function:

QwTrack::QwTrack ( const QwTrack that)

Copy constructor from pointer.

Copy constructor from pointer

Parameters
thatOriginal track

Definition at line 59 of file QwTrack.cc.

References fBack, fFront, and Initialize().

60 : VQwTrackingElement(*that)
61 {
62  // Initialize
63  Initialize();
64 
65  // Null pointer
66  if (that == 0) return;
67 
68  *this = *that;
69 
70  // Copy tracks
71  if (that->fFront) fFront = new QwPartialTrack(that->fFront);
72  if (that->fBack) fBack = new QwPartialTrack(that->fBack);
73 }
void Initialize()
Initialization.
Definition: QwTrack.cc:149
VQwTrackingElement()
Default constructor.
QwPartialTrack * fBack
Back partial track.
Definition: QwTrack.h:148
QwPartialTrack * fFront
Front partial track.
Definition: QwTrack.h:147
Contains the straight part of a track in one region only.

+ Here is the call graph for this function:

QwTrack::~QwTrack ( )
virtual

Virtual destructor.

Virtual destructor

Definition at line 78 of file QwTrack.cc.

References ClearPartialTracks(), fBack, and fFront.

79 {
80  // Delete objects
81  if (fFront) delete fFront;
82  if (fBack) delete fBack;
83 
85 }
QwPartialTrack * fBack
Back partial track.
Definition: QwTrack.h:148
QwPartialTrack * fFront
Front partial track.
Definition: QwTrack.h:147
void ClearPartialTracks(Option_t *option="")
Definition: QwTrack.cc:186

+ Here is the call graph for this function:

Member Function Documentation

void QwTrack::AddPartialTrack ( const QwPartialTrack partialtrack)

Definition at line 171 of file QwTrack.cc.

References fNQwPartialTracks, and fQwPartialTracks.

Referenced by AddPartialTrackList(), CreateNewPartialTrack(), and QwTrack().

172 {
173  fQwPartialTracks.push_back(new QwPartialTrack(partialtrack));
175 }
Int_t fNQwPartialTracks
Number of partial tracks in this track.
Definition: QwTrack.h:35
std::vector< QwPartialTrack * > fQwPartialTracks
List of partial tracks in this track.
Definition: QwTrack.h:37
Contains the straight part of a track in one region only.

+ Here is the caller graph for this function:

void QwTrack::AddPartialTrackList ( const std::vector< QwPartialTrack * > &  partialtracklist)

Definition at line 178 of file QwTrack.cc.

References AddPartialTrack().

Referenced by operator=().

179 {
180  for (std::vector<QwPartialTrack*>::const_iterator partialtrack = partialtracklist.begin();
181  partialtrack != partialtracklist.end(); partialtrack++)
182  AddPartialTrack(*partialtrack);
183 }
void AddPartialTrack(const QwPartialTrack *partialtrack)
Definition: QwTrack.cc:171

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

QwTrack::ClassDef ( QwTrack  ,
 
)
void QwTrack::ClearPartialTracks ( Option_t *  option = "")

Definition at line 186 of file QwTrack.cc.

References fNQwPartialTracks, and fQwPartialTracks.

Referenced by operator=(), ResetPartialTracks(), and ~QwTrack().

187 {
188  for (size_t i = 0; i < fQwPartialTracks.size(); i++) {
189  QwPartialTrack* tl = fQwPartialTracks.at(i);
190  delete tl;
191  }
192  fQwPartialTracks.clear();
193  fNQwPartialTracks = 0;
194 }
Int_t fNQwPartialTracks
Number of partial tracks in this track.
Definition: QwTrack.h:35
std::vector< QwPartialTrack * > fQwPartialTracks
List of partial tracks in this track.
Definition: QwTrack.h:37
Contains the straight part of a track in one region only.

+ Here is the caller graph for this function:

QwPartialTrack * QwTrack::CreateNewPartialTrack ( )

Creating and adding partial tracks

Definition at line 163 of file QwTrack.cc.

References AddPartialTrack().

164 {
165  QwPartialTrack* partialtrack = new QwPartialTrack();
166  AddPartialTrack(partialtrack);
167  return partialtrack;
168 }
void AddPartialTrack(const QwPartialTrack *partialtrack)
Definition: QwTrack.cc:171
Contains the straight part of a track in one region only.

+ Here is the call graph for this function:

const TVector3& QwTrack::GetMagneticFieldIntegral ( ) const
inline

Get magnetic field integral.

Definition at line 81 of file QwTrack.h.

References fBdlXYZ.

81  {
82  return fBdlXYZ;
83  }
TVector3 fBdlXYZ
Magnetic field integral along track.
Definition: QwTrack.h:144
Int_t QwTrack::GetNumberOfPartialTracks ( ) const
inline

Definition at line 67 of file QwTrack.h.

References fNQwPartialTracks.

67 { return fNQwPartialTracks; };
Int_t fNQwPartialTracks
Number of partial tracks in this track.
Definition: QwTrack.h:35
const QwPartialTrack* QwTrack::GetPartialTrack ( const int  pt) const
inline

Get the specified partial track.

Definition at line 69 of file QwTrack.h.

References fQwPartialTracks.

69 { return fQwPartialTracks.at(pt); };
std::vector< QwPartialTrack * > fQwPartialTracks
List of partial tracks in this track.
Definition: QwTrack.h:37
void QwTrack::Initialize ( )

Initialization.

Initialization

Definition at line 149 of file QwTrack.cc.

References fBack, fChi, fFront, and fMomentum.

Referenced by QwTrack().

150 {
151  // Initialize the members;
152  fChi = 0.0;
153  fMomentum = 0.0;
154 
155  // Initialize all pointers
156  fFront = 0;
157  fBack = 0;
158 }
double fChi
Combined chi, i.e. the sum of the chi for front and back partial track.
Definition: QwTrack.h:99
QwPartialTrack * fBack
Back partial track.
Definition: QwTrack.h:148
QwPartialTrack * fFront
Front partial track.
Definition: QwTrack.h:147
double fMomentum
Spectrometer momentum.
Definition: QwTrack.h:100

+ Here is the caller graph for this function:

QwTrack & QwTrack::operator= ( const QwTrack that)

Assignment operator.

Assignment operator

Parameters
thatOriginal track

Definition at line 91 of file QwTrack.cc.

References AddPartialTrackList(), ClearPartialTracks(), fChi, fDirectionDiff, fDirectionPhioff, fDirectionThetaoff, fDirectionXoff, fDirectionYoff, fDirectionZoff, fEndDirectionActual, fEndDirectionActualRK4, fEndDirectionActualRKF45, fEndDirectionGoal, fEndPositionActual, fEndPositionActualRK4, fEndPositionActualRKF45, fEndPositionGoal, fIterationsNewton, fIterationsRK4, fIterationsRKF45, fIterationsRungeKutta, fMomentum, fPhi, fPositionDiff, fPositionPhioff, fPositionRoff, fPositionThetaoff, fPositionXoff, fPositionYoff, fQwPartialTracks, fStartDirection, fStartPosition, fTheta, fVertexR, fVertexZ, and VQwTrackingElement::operator=().

92 {
93  if (this == &that) return *this;
94 
96 
97  fPhi = that.fPhi;
98  fTheta = that.fTheta;
99  fVertexZ = that.fVertexZ;
100  fVertexR = that.fVertexR;
101 
102  fChi = that.fChi;
103  fMomentum = that.fMomentum;
104 
107 
114 
121 
124 
127 
130 
134 
138 
139  // Copy partial tracks
142 
143  return *this;
144 }
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
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
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 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
double fDirectionPhioff
Difference in momentum azimuthal angle at matching plane.
Definition: QwTrack.h:117
double fPositionRoff
Difference in radial position at matching plane.
Definition: QwTrack.h:110
TVector3 fEndDirectionActual
Actual direction of track at back plane.
Definition: QwTrack.h:130
TVector3 fStartDirection
Start direction of front track.
Definition: QwTrack.h:124
VQwTrackingElement & operator=(const VQwTrackingElement &that)
Assignment operator.
TVector3 fEndPositionActual
Actual position of track at back plane.
Definition: QwTrack.h:129
double fDirectionXoff
Difference in X momentum at matching plane.
Definition: QwTrack.h:114
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
TVector3 fPositionDiff
Matching of front and back track position and direction at matching plane.
Definition: QwTrack.h:107
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
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
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
double fVertexZ
Primary vertex position in longitudinal direction.
Definition: QwTrack.h:95
void ClearPartialTracks(Option_t *option="")
Definition: QwTrack.cc:186
double fVertexR
Primary vertex position in transverse direction.
Definition: QwTrack.h:96

+ Here is the call graph for this function:

void QwTrack::PrintPartialTracks ( Option_t *  option = "") const

Definition at line 203 of file QwTrack.cc.

References QwLog::endl(), fQwPartialTracks, and QwMessage.

204 {
205  for (std::vector<QwPartialTrack*>::const_iterator partialtrack = fQwPartialTracks.begin();
206  partialtrack != fQwPartialTracks.end(); partialtrack++) {
207  QwMessage << **partialtrack << QwLog::endl;
208  }
209 }
#define QwMessage
Predefined log drain for regular messages.
Definition: QwLog.h:50
static std::ostream & endl(std::ostream &)
End of the line.
Definition: QwLog.cc:299
std::vector< QwPartialTrack * > fQwPartialTracks
List of partial tracks in this track.
Definition: QwTrack.h:37

+ Here is the call graph for this function:

void QwTrack::ResetPartialTracks ( Option_t *  option = "")

Definition at line 197 of file QwTrack.cc.

References ClearPartialTracks().

198 {
200 }
void ClearPartialTracks(Option_t *option="")
Definition: QwTrack.cc:186

+ Here is the call graph for this function:

void QwTrack::SetMagneticFieldIntegral ( double  bdl,
double  bdlx,
double  bdly,
double  bdlz 
)
inline

Set magnetic field integral.

Definition at line 76 of file QwTrack.h.

References fBdl, and fBdlXYZ.

Referenced by QwRayTracer::Bridge().

76  {
77  fBdl = bdl;
78  fBdlXYZ = TVector3(bdlx,bdly,bdlz);
79  }
double fBdl
Magnetic field integrals.
Definition: QwTrack.h:143
TVector3 fBdlXYZ
Magnetic field integral along track.
Definition: QwTrack.h:144

+ Here is the caller graph for this function:

Friends And Related Function Documentation

std::ostream& operator<< ( std::ostream &  stream,
const QwTrack t 
)
friend

Output stream operator for tracks.

Output stream operator overloading

Definition at line 215 of file QwTrack.cc.

216 {
217  stream << "track: ";
218  stream << "(" << t.GetRegion() << "/" << "?UD"[t.GetPackage()] << ") ";
219  stream << "Start: " << t.fStartPosition << "/" << t.fStartDirection << std::endl;
220  stream << "End (goal): " << t.fEndPositionGoal << "/" << t.fEndDirectionGoal << std::endl;
221  stream << "End (actual): " << t.fEndPositionActual << "/" << t.fEndDirectionActual << std::endl;
222  return stream;
223 }
TVector3 fStartPosition
Position and direction before and after swimming.
Definition: QwTrack.h:123
EQwRegionID GetRegion() const
Get the region.
TVector3 fEndDirectionActual
Actual direction of track at back plane.
Definition: QwTrack.h:130
TVector3 fStartDirection
Start direction of front track.
Definition: QwTrack.h:124
TVector3 fEndPositionActual
Actual position of track at back plane.
Definition: QwTrack.h:129
TVector3 fEndDirectionGoal
Goal direction of back track.
Definition: QwTrack.h:127
TVector3 fEndPositionGoal
Goal position of back track.
Definition: QwTrack.h:126
EQwDetectorPackage GetPackage() const
Get the package.

Field Documentation

QwPartialTrack* QwTrack::fBack

Back partial track.

Definition at line 148 of file QwTrack.h.

Referenced by DrawTreeLineVDC(), Initialize(), QwTrack(), and ~QwTrack().

double QwTrack::fBdl

Magnetic field integrals.

Magnetic field integral along track

Definition at line 143 of file QwTrack.h.

Referenced by SetMagneticFieldIntegral().

TVector3 QwTrack::fBdlXYZ

Magnetic field integral along track.

Definition at line 144 of file QwTrack.h.

Referenced by GetMagneticFieldIntegral(), and SetMagneticFieldIntegral().

double QwTrack::fChi

Combined chi, i.e. the sum of the chi for front and back partial track.

Definition at line 99 of file QwTrack.h.

Referenced by QwRayTracer::Bridge(), Initialize(), and operator=().

TVector3 QwTrack::fDirectionDiff

Difference in momentum at matching plane.

Definition at line 113 of file QwTrack.h.

Referenced by QwRayTracer::Bridge(), and operator=().

double QwTrack::fDirectionPhioff

Difference in momentum azimuthal angle at matching plane.

Definition at line 117 of file QwTrack.h.

Referenced by QwRayTracer::Bridge(), and operator=().

double QwTrack::fDirectionThetaoff

Difference in momentum polar angle at matching plane.

Definition at line 118 of file QwTrack.h.

Referenced by QwRayTracer::Bridge(), and operator=().

double QwTrack::fDirectionXoff

Difference in X momentum at matching plane.

Definition at line 114 of file QwTrack.h.

Referenced by QwRayTracer::Bridge(), and operator=().

double QwTrack::fDirectionYoff

Difference in Y momentum at matching plane.

Definition at line 115 of file QwTrack.h.

Referenced by QwRayTracer::Bridge(), and operator=().

double QwTrack::fDirectionZoff

Difference in Z momentum at matching plane.

Definition at line 116 of file QwTrack.h.

Referenced by QwRayTracer::Bridge(), and operator=().

TVector3 QwTrack::fEndDirectionActual

Actual direction of track at back plane.

Definition at line 130 of file QwTrack.h.

Referenced by QwMatrixLookup::Bridge(), QwRayTracer::Bridge(), operator<<(), and operator=().

TVector3 QwTrack::fEndDirectionActualRK4

Actual direction of track at back plane using Runge-Kutta 4th order.

Definition at line 134 of file QwTrack.h.

Referenced by QwRayTracer::Bridge(), and operator=().

TVector3 QwTrack::fEndDirectionActualRKF45

Actual direction of track at back plane using Runge-Kutta-Fehlberg.

Definition at line 138 of file QwTrack.h.

Referenced by QwRayTracer::Bridge(), and operator=().

TVector3 QwTrack::fEndDirectionGoal

Goal direction of back track.

Definition at line 127 of file QwTrack.h.

Referenced by QwMatrixLookup::Bridge(), QwRayTracer::Bridge(), operator<<(), and operator=().

TVector3 QwTrack::fEndPositionActual

Actual position of track at back plane.

Definition at line 129 of file QwTrack.h.

Referenced by QwMatrixLookup::Bridge(), QwRayTracer::Bridge(), operator<<(), and operator=().

TVector3 QwTrack::fEndPositionActualRK4

Actual position of track at back plane using Runge-Kutta 4th order.

Definition at line 133 of file QwTrack.h.

Referenced by QwRayTracer::Bridge(), and operator=().

TVector3 QwTrack::fEndPositionActualRKF45

Actual position of track at back plane using Runge-Kutta-Fehlberg.

Definition at line 137 of file QwTrack.h.

Referenced by QwRayTracer::Bridge(), and operator=().

TVector3 QwTrack::fEndPositionGoal

Goal position of back track.

Definition at line 126 of file QwTrack.h.

Referenced by QwMatrixLookup::Bridge(), QwRayTracer::Bridge(), operator<<(), and operator=().

QwPartialTrack* QwTrack::fFront

Front partial track.

Definition at line 147 of file QwTrack.h.

Referenced by QwEvent::CalculateKinematics(), Initialize(), QwTrack(), and ~QwTrack().

int QwTrack::fIterationsNewton

Number of iterations in Newton's method.

Definition at line 102 of file QwTrack.h.

Referenced by QwRayTracer::Bridge(), and operator=().

int QwTrack::fIterationsRK4

Number of iterations using Runge-Kutta 4th order.

Definition at line 132 of file QwTrack.h.

Referenced by QwRayTracer::Bridge(), and operator=().

int QwTrack::fIterationsRKF45

Number of iterations using Runge-Kutta-Fehlberg.

Definition at line 136 of file QwTrack.h.

Referenced by QwRayTracer::Bridge(), and operator=().

int QwTrack::fIterationsRungeKutta

Number of iterations in Runge-Kutta method.

Definition at line 103 of file QwTrack.h.

Referenced by QwRayTracer::Bridge(), and operator=().

double QwTrack::fMomentum

Spectrometer momentum.

Definition at line 100 of file QwTrack.h.

Referenced by QwMatrixLookup::Bridge(), QwRayTracer::Bridge(), QwEvent::CalculateKinematics(), Initialize(), and operator=().

Int_t QwTrack::fNQwPartialTracks
private

Number of partial tracks in this track.

Definition at line 35 of file QwTrack.h.

Referenced by AddPartialTrack(), ClearPartialTracks(), and GetNumberOfPartialTracks().

double QwTrack::fPhi

Quantities determined from front partial track.

Azimuthal angle phi of track at primary vertex

Definition at line 93 of file QwTrack.h.

Referenced by operator=().

TVector3 QwTrack::fPositionDiff

Matching of front and back track position and direction at matching plane.

Difference in position at matching plane

Definition at line 107 of file QwTrack.h.

Referenced by QwRayTracer::Bridge(), and operator=().

double QwTrack::fPositionPhioff

Difference in azimuthal angle phi at matching plane.

Definition at line 111 of file QwTrack.h.

Referenced by QwRayTracer::Bridge(), and operator=().

double QwTrack::fPositionRoff

Difference in radial position at matching plane.

Definition at line 110 of file QwTrack.h.

Referenced by QwRayTracer::Bridge(), and operator=().

double QwTrack::fPositionThetaoff

Difference in polar angle theta at matching plane.

Definition at line 112 of file QwTrack.h.

Referenced by QwRayTracer::Bridge(), and operator=().

double QwTrack::fPositionXoff

Difference in X position at matching plane.

Definition at line 108 of file QwTrack.h.

Referenced by QwRayTracer::Bridge(), and operator=().

double QwTrack::fPositionYoff

Difference in Y position at matching plane.

Definition at line 109 of file QwTrack.h.

Referenced by QwRayTracer::Bridge(), and operator=().

std::vector<QwPartialTrack*> QwTrack::fQwPartialTracks
private

List of partial tracks in this track.

Definition at line 37 of file QwTrack.h.

Referenced by AddPartialTrack(), ClearPartialTracks(), GetPartialTrack(), operator=(), and PrintPartialTracks().

TVector3 QwTrack::fStartDirection

Start direction of front track.

Definition at line 124 of file QwTrack.h.

Referenced by QwMatrixLookup::Bridge(), QwRayTracer::Bridge(), operator<<(), and operator=().

TVector3 QwTrack::fStartPosition

Position and direction before and after swimming.

Start position of front track

Definition at line 123 of file QwTrack.h.

Referenced by QwMatrixLookup::Bridge(), QwRayTracer::Bridge(), operator<<(), and operator=().

double QwTrack::fTheta

Polar angle theta of track at primary vertex.

Definition at line 94 of file QwTrack.h.

Referenced by operator=().

double QwTrack::fVertexR

Primary vertex position in transverse direction.

Definition at line 96 of file QwTrack.h.

Referenced by operator=().

double QwTrack::fVertexZ

Primary vertex position in longitudinal direction.

Definition at line 95 of file QwTrack.h.

Referenced by operator=().


The documentation for this class was generated from the following files: