FEDRA emulsion software from the OPERA Collaboration
VERTEX::Track Class Reference

#include <VtTrack.hh>

Inheritance diagram for VERTEX::Track:
Collaboration diagram for VERTEX::Track:

Public Member Functions

— Constructors —
 Track ()
 
 Track (const Track *const track)
 
 Track (const Track &rhs)
 
 Track (const MATRIX::VtVector &v, const MATRIX::CMatrix &c)
 
virtual ~Track ()
 
void set (double x, double y, double z, double tx, double ty, double p, const MATRIX::CMatrix &c)
 
Trackoperator= (const Track &rhs)
 
— TrackIf access methods —
float x () const
 x at z=0 (in Vt: p(3,..)) More...
 
float y () const
 y at z=0 (in Vt: p(4,..)) More...
 
float x (float z) const
 Track $x$ position at $z$. More...
 
float y (float z) const
 Track $y$ position at $z$. More...
 
float z () const
 z = 0 in Vt More...
 
float tx () const
 slope (in Vt: p(1,..)) More...
 
float ty () const
 slope (in Vt: p(2,..)) More...
 
float p () const
 momentum (in Vt: p(5,..)) More...
 
float pt () const
 transv. momentum $\sqrt{p_x^2 + p_y^2}$ More...
 
float pz () const
 $p \cdot e_z$ More...
 
float chi2 () const
 dummy function: always return 0 More...
 
unsigned short int ndf () const
 dummy function: always return 0 More...
 
float phi () const
 azimuthal angle $\phi$ More...
 
float theta () const
 polar angle $\theta = \cos^{-1}(e_z)$ More...
 
float eta () const
 rapidity $-\log\tan(\theta/2.)$ More...
 
int charge () const
 charge More...
 
float energy (double rm=0.) const
 Energy with given rest-mass (in GeV) $E = \sqrt{p^2 + m^2}$. More...
 
float xf (double rm=0.) const
 Feynman Variable $x_F = \frac{E+p_z}{(E+p_z)_{max}}$. More...
 
float rap (double rm=0.) const
 Rapidity $y = \frac{1}{2}\ln(\frac{E+p_z}{E-p_z})$. More...
 
MATRIX::VtVector evec () const
 $\vec{v} = (e_x,e_y,e_z)$ unit vector along refitted track More...
 
MATRIX::VtVector tvec () const
 $\vec{v} = (t_x,t_y,1.)$ refitted slope vector More...
 
MATRIX::VtVector pvec () const
 $\vec{v} = (p_x,p_y,p_z)$ refitted mom. vector More...
 
float cov_x (double dz=0.) const
 get $\sigma_x^2|_{z+dz}$ More...
 
float cov_y (double dz=0.) const
 get $\sigma_y^2|_{z+dz}$ More...
 
float cov_tx () const
 get $\sigma_{t_x}^2$ More...
 
float cov_ty () const
 get $\sigma_{t_y}^2$ More...
 
float cov_p () const
 get $\sigma_p^2$ More...
 
bool isValid () const
 dummy function: always return true More...
 
void valid ()
 dummy function: do nothing More...
 
void invalid ()
 dummy function: do nothing More...
 
bool propagate (const double zz)
 propagate track to $z$ More...
 
— Other access methods —
double px () const
 $p \cdot e_x$ More...
 
double py () const
 $p \cdot e_y$ More...
 
double ex () const
 $t_x \cdot e_z$ More...
 
double ey () const
 $t_y \cdot e_z$ More...
 
double ez () const
 $1/\sqrt{1+t_x^2+t_y^2}$ More...
 
double xerr (double dz=0) const
 $\sqrt{\sigma_x^2(z+dz)}$ More...
 
double yerr (double dz=0) const
 $\sqrt{\sigma_y^2(z+dz)}$ More...
 
double txerr () const
 $\sqrt{\sigma^2_{tx}}$ More...
 
double tyerr () const
 $\sqrt{\sigma^2_{ty}}$ More...
 
double perr () const
 $\sqrt{\sigma^2_p}$ More...
 
const MATRIX::CMatrixV () const
 covariance matrix More...
 
const MATRIX::CMatrixG () const
 inverse $5\times5$ covariance matrix More...
 
const MATRIX::CMatrixGM () const
 inverse $4\times4$ covariance matrix, without momentum More...
 
— Expert functions —
std::ostream & print (std::ostream &os) const
 called by cout More...
 
double rm () const
 rest-mass More...
 
void rm (const double mass)
 set rest-mass (needed for mass constrained fits) More...
 
double rmCC () const
 conjugated rest-mass More...
 
void rmCC (const double mass)
 set conjugated rest-mass More...
 
void delete_mom ()
 
- Public Member Functions inherited from VERTEX::RelationList
 RelationList ()
 
virtual ~RelationList ()
 
virtual void clear ()
 clear all relations More...
 
void push__back (Relation *rel)
 add a relation More...
 
Relationback () const
 return last relation More...
 
unsigned int size () const
 return no of relations More...
 
virtual const iterator erase (const iterator &pos)
 Erase an element. More...
 
iterator begin ()
 
iterator end ()
 
const_iterator begin () const
 
const_iterator end () const
 
reverse_iterator rbegin ()
 
reverse_iterator rend ()
 
const_reverse_iterator rbegin () const
 
const_reverse_iterator rend () const
 
Relation *const operator[] (unsigned int i) const
 direct access More...
 
void print (std::ostream &os) const
 
void unbook (Relation *const rel)
 get rid of relation pointer More...
 
virtual void remove (Relation *const rel)
 remove a relation More...
 
const RelationListoperator= (const RelationList &rhs)
 
bool operator== (const RelationList &rhs) const
 

— Operators —

short int t_Q
 
double t_rm
 
double t_rmCC
 
MATRIX::VtVector t_p
 
MATRIX::CMatrix t_V
 
MATRIX::CMatrix t_G
 
MATRIX::CMatrix t_GM
 
 ClassDef (Track, 1)
 

Additional Inherited Members

- Public Types inherited from VERTEX::RelationList
typedef ConstRelationIterator const_iterator
 
typedef ConstReverseRelationIterator const_reverse_iterator
 
typedef RelationIterator iterator
 
typedef ReverseRelationIterator reverse_iterator
 

Detailed Description

Track representation class

Constructor & Destructor Documentation

◆ Track() [1/4]

Track::Track ( )
64 :
65 t_Q (0),
66 t_rm (0.),
67 t_rmCC (0.),
68 t_p (VtVector(6))
69 {}
Definition: VtVector.hh:45
double t_rmCC
Definition: VtTrack.hh:226
double t_rm
Definition: VtTrack.hh:225
MATRIX::VtVector t_p
Definition: VtTrack.hh:228
short int t_Q
Definition: VtTrack.hh:224

◆ Track() [2/4]

Track::Track ( const Track *const  track)
71 :
72 t_Q (track->t_Q),
73 t_rm (track->t_rm),
74 t_rmCC (track->t_rmCC),
75 t_p (track->t_p),
76 t_V (track->t_V)
77 {
78 propagate(0.); // propagate track to z=0.
79 t_G = t_V.dsinv();
80 t_GM.copy(t_V.dsinv(4)); // inverse without momentum
81 }
void copy(const VtSymMatrix &rhs)
in case matrix dimensions differ
Definition: VtSymMatrix.C:649
const VtSymMatrix dsinv(int dim=0) const
return inverse
Definition: VtSymMatrix.C:275
MATRIX::CMatrix t_V
Definition: VtTrack.hh:230
bool propagate(const double zz)
propagate track to $z$
Definition: VtTrack.C:230
MATRIX::CMatrix t_G
Definition: VtTrack.hh:231
MATRIX::CMatrix t_GM
Definition: VtTrack.hh:232
Definition: bitview.h:14

◆ Track() [3/4]

Track::Track ( const Track rhs)
83 :
84 t_Q (rhs.t_Q),
85 t_rm (rhs.t_rm),
86 t_rmCC (rhs.t_rmCC),
87 t_p (rhs.t_p),
88 t_V (rhs.t_V)
89 {
90 propagate(0.); // propagate track to z=0.
91 t_G = t_V.dsinv();
92 t_GM.copy(t_V.dsinv(4)); // inverse without momentum
93 }

◆ Track() [4/4]

VERTEX::Track::Track ( const MATRIX::VtVector v,
const MATRIX::CMatrix c 
)

◆ ~Track()

Track::~Track ( )
virtual
95{}

Member Function Documentation

◆ charge()

int Track::charge ( ) const

charge

163{ return t_Q; }

◆ chi2()

float Track::chi2 ( ) const

dummy function: always return 0

158{ return 0; }

◆ ClassDef()

VERTEX::Track::ClassDef ( Track  ,
 
)
protected

◆ cov_p()

float Track::cov_p ( ) const

get $\sigma_p^2$

194{ return t_V.p(); }
double p() const
rtra->cvf[14] $* p^4$
Definition: CMatrix.hh:111

◆ cov_tx()

float Track::cov_tx ( ) const

get $\sigma_{t_x}^2$

192{ return t_V(2,2); }

◆ cov_ty()

float Track::cov_ty ( ) const

get $\sigma_{t_y}^2$

193{ return t_V(3,3); }

◆ cov_x()

float Track::cov_x ( double  dz = 0.) const

get $\sigma_x^2|_{z+dz}$

184 {
185 return t_V.x_prop(dz);
186 }
brick dz
Definition: RecDispMC.C:107
double x_prop(double dz) const
Definition: CMatrix.hh:147

◆ cov_y()

float Track::cov_y ( double  dz = 0.) const

get $\sigma_y^2|_{z+dz}$

188 {
189 return t_V.y_prop(dz);
190 }
double y_prop(double dz) const
Definition: CMatrix.hh:157

◆ delete_mom()

void Track::delete_mom ( )
244 {
245 t_V.set_x_p(0.);
246 t_V.set_y_p(0.);
247 t_V.set_tx_p(0.);
248 t_V.set_ty_p(0.);
249 t_V.set_p(1e5);
250 t_G = t_V;
251 if(t_G.VtDsinv()==false)
252 cout << " error inverting matrix!!" << endl;
253
254 }
void set_ty_p(const double typ)
Definition: CMatrix.hh:141
void set_tx_p(const double txp)
Definition: CMatrix.hh:137
void set_x_p(const double xp)
Definition: CMatrix.hh:123
void set_p(const double p)
Definition: CMatrix.hh:143
void set_y_p(const double yp)
Definition: CMatrix.hh:131
bool VtDsinv(int dim=0)
transform to inverse
Definition: VtSymMatrix.C:292

◆ energy()

float Track::energy ( double  rm = 0.) const

Energy with given rest-mass (in GeV) $E = \sqrt{p^2 + m^2}$.

164{ return sqrt(sqr(t_p[5]) + sqr(rm)); }
double rm() const
rest-mass
const T sqr(const T &x)
compute the square of a number: $x*x$
Definition: VtUtil.hh:37

◆ eta()

float Track::eta ( ) const

rapidity $-\log\tan(\theta/2.)$

162{ return -log(tan(theta()/2.)); }
FILE * log
float theta() const
polar angle $\theta = \cos^{-1}(e_z)$
Definition: VtTrack.C:161

◆ evec()

MATRIX::VtVector VERTEX::Track::evec ( ) const
inline

$\vec{v} = (e_x,e_y,e_z)$ unit vector along refitted track

137{ return MATRIX::VtVector(ex(),ey(),ez()); }
double ez() const
$1/\sqrt{1+t_x^2+t_y^2}$
double ey() const
$t_y \cdot e_z$
double ex() const
$t_x \cdot e_z$

◆ ex()

double VERTEX::Track::ex ( ) const

$t_x \cdot e_z$

◆ ey()

double VERTEX::Track::ey ( ) const

$t_y \cdot e_z$

◆ ez()

double VERTEX::Track::ez ( ) const

$1/\sqrt{1+t_x^2+t_y^2}$

◆ G()

const MATRIX::CMatrix & Track::G ( ) const

inverse $5\times5$ covariance matrix

206{ return t_G; }

◆ GM()

const MATRIX::CMatrix & Track::GM ( ) const

inverse $4\times4$ covariance matrix, without momentum

207{ return t_GM; }

◆ invalid()

void Track::invalid ( )

dummy function: do nothing

197{ return; }

◆ isValid()

bool Track::isValid ( ) const

dummy function: always return true

195{ return true; }

◆ ndf()

unsigned short int Track::ndf ( ) const

dummy function: always return 0

159{ return 0; }

◆ operator=()

Track & Track::operator= ( const Track rhs)
128 {
129 if(this == &rhs) { return *this; }
130
132
133 t_Q = rhs.t_Q;
134 t_rm = rhs.t_rm;
135 t_rmCC = rhs.t_rmCC;
136
137 t_p = rhs.t_p;
138 t_V = rhs.t_V;
139 t_G = rhs.t_G;
140 t_GM = rhs.t_GM;
141
142 return *this;
143 }
const RelationList & operator=(const RelationList &rhs)
Definition: VtRelationList.C:121

◆ p()

float Track::p ( ) const

momentum (in Vt: p(5,..))

155{ return t_p[5]; }

◆ perr()

double VERTEX::Track::perr ( ) const

$\sqrt{\sigma^2_p}$

◆ phi()

float Track::phi ( ) const

azimuthal angle $\phi$

160{ return (180. + atan2(t_p[3],t_p[4])*57.29577957); }

◆ print()

std::ostream & Track::print ( std::ostream &  os) const

called by cout

213 {
214
215 return os << " x: " << x()
216 << " y: " << y()
217 << " z: " << z()
218 << " tx: " << tx()
219 << " ty: " << ty()
220 << " p: " << p()
221 << endl
222 << " Cov Matrix: " << V()
223 << " inverse: " << G()
224 << " inverse without momentum: " << GM();
225 }
float x() const
x at z=0 (in Vt: p(3,..))
Definition: VtTrack.C:148
float z() const
z = 0 in Vt
Definition: VtTrack.C:150
const MATRIX::CMatrix & GM() const
inverse $4\times4$ covariance matrix, without momentum
Definition: VtTrack.C:207
float ty() const
slope (in Vt: p(2,..))
Definition: VtTrack.C:154
float y() const
y at z=0 (in Vt: p(4,..))
Definition: VtTrack.C:149
const MATRIX::CMatrix & G() const
inverse $5\times5$ covariance matrix
Definition: VtTrack.C:206
const MATRIX::CMatrix & V() const
covariance matrix
Definition: VtTrack.C:205
float p() const
momentum (in Vt: p(5,..))
Definition: VtTrack.C:155
float tx() const
slope (in Vt: p(1,..))
Definition: VtTrack.C:153

◆ propagate()

bool Track::propagate ( const double  zz)

propagate track to $z$

230 {
231 if( zz == t_p[2] ) return true;
232
233 const double dz = zz - t_p[2];
234 t_p[0] += dz * t_p[3];
235 t_p[1] += dz * t_p[4];
236 t_p[2] = zz;
237 t_V.propagate(dz); // propagate Track covariance matrix
238 return true;
239 }
void propagate(const double dz)
propagate all matrix elements
Definition: CMatrix.C:74

◆ pt()

float Track::pt ( ) const

transv. momentum $\sqrt{p_x^2 + p_y^2}$

156{ return sqrt(sqr(px()) + sqr(py())); }
double py() const
$p \cdot e_y$
double px() const
$p \cdot e_x$

◆ pvec()

MATRIX::VtVector VERTEX::Track::pvec ( ) const
inline

$\vec{v} = (p_x,p_y,p_z)$ refitted mom. vector

141{ return MATRIX::VtVector(px(),py(),pz()); }
float pz() const
$p \cdot e_z$
Definition: VtTrack.C:157

◆ px()

double VERTEX::Track::px ( ) const

$p \cdot e_x$

◆ py()

double VERTEX::Track::py ( ) const

$p \cdot e_y$

◆ pz()

float Track::pz ( ) const

$p \cdot e_z$

157{ return t_p[5] * ez(); }

◆ rap()

float Track::rap ( double  rm = 0.) const

Rapidity $y = \frac{1}{2}\ln(\frac{E+p_z}{E-p_z})$.

166{ return 0.5*log((energy(mass)+pz())/(energy(mass)-pz())); }
float energy(double rm=0.) const
Energy with given rest-mass (in GeV) $E = \sqrt{p^2 + m^2}$.
Definition: VtTrack.C:164
float mass
Definition: check_vertex.C:21

◆ rm() [1/2]

double VERTEX::Track::rm ( ) const

rest-mass

◆ rm() [2/2]

void VERTEX::Track::rm ( const double  mass)

set rest-mass (needed for mass constrained fits)

◆ rmCC() [1/2]

double VERTEX::Track::rmCC ( ) const

conjugated rest-mass

◆ rmCC() [2/2]

void VERTEX::Track::rmCC ( const double  mass)

set conjugated rest-mass

◆ set()

void Track::set ( double  x,
double  y,
double  z,
double  tx,
double  ty,
double  p,
const MATRIX::CMatrix c 
)
112 {
113 t_p[0] = x;
114 t_p[1] = y;
115 t_p[2] = z;
116 t_p[3] = tx;
117 t_p[4] = ty;
118 t_p[5] = p;
119 t_V.copy(c);
120 propagate(0.); // propagate track to z=0.
121 t_G = t_V.dsinv();
122 t_GM.copy(t_V.dsinv(4)); // inverse without momentum
123 }

◆ theta()

float Track::theta ( ) const

polar angle $\theta = \cos^{-1}(e_z)$

161{ return acos(ez()); }

◆ tvec()

MATRIX::VtVector VERTEX::Track::tvec ( ) const
inline

$\vec{v} = (t_x,t_y,1.)$ refitted slope vector

139{ return MATRIX::VtVector(tx(),ty(),1.); }

◆ tx()

float Track::tx ( ) const

slope (in Vt: p(1,..))

153{ return t_p[3]; }

◆ txerr()

double VERTEX::Track::txerr ( ) const

$\sqrt{\sigma^2_{tx}}$

◆ ty()

float Track::ty ( ) const

slope (in Vt: p(2,..))

154{ return t_p[4]; }

◆ tyerr()

double VERTEX::Track::tyerr ( ) const

$\sqrt{\sigma^2_{ty}}$

◆ V()

const MATRIX::CMatrix & Track::V ( ) const

covariance matrix

205{ return t_V; }

◆ valid()

void Track::valid ( )

dummy function: do nothing

196{ return; }

◆ x() [1/2]

float Track::x ( ) const

x at z=0 (in Vt: p(3,..))

148{ return t_p[0]; }

◆ x() [2/2]

float Track::x ( float  z) const

Track $x$ position at $z$.

151{ return t_p[0] + t_p[3] * (z-t_p[2]); }

◆ xerr()

double VERTEX::Track::xerr ( double  dz = 0) const

$\sqrt{\sigma_x^2(z+dz)}$

◆ xf()

float Track::xf ( double  rm = 0.) const

Feynman Variable $x_F = \frac{E+p_z}{(E+p_z)_{max}}$.

165{ return (pz() - energy(mass)*0.99898)/0.93827; }

◆ y() [1/2]

float Track::y ( ) const

y at z=0 (in Vt: p(4,..))

149{ return t_p[1]; }

◆ y() [2/2]

float Track::y ( float  z) const

Track $y$ position at $z$.

152{ return t_p[1] + t_p[4] * (z-t_p[2]); }

◆ yerr()

double VERTEX::Track::yerr ( double  dz = 0) const

$\sqrt{\sigma_y^2(z+dz)}$

◆ z()

float Track::z ( ) const

z = 0 in Vt

150{ return t_p[2]; }

Member Data Documentation

◆ t_G

MATRIX::CMatrix VERTEX::Track::t_G
protected

◆ t_GM

MATRIX::CMatrix VERTEX::Track::t_GM
protected

◆ t_p

MATRIX::VtVector VERTEX::Track::t_p
protected

◆ t_Q

short int VERTEX::Track::t_Q
private

◆ t_rm

double VERTEX::Track::t_rm
private

◆ t_rmCC

double VERTEX::Track::t_rmCC
private

◆ t_V

MATRIX::CMatrix VERTEX::Track::t_V
protected

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