FEDRA emulsion software from the OPERA Collaboration
VtTrack.hh
Go to the documentation of this file.
1#ifndef __VTTRACK_HH
2#define __VTTRACK_HH
3// *****************************************************************************
4//
5// source:
6//
7// type: source code
8//
9// created: 21. Aug 2000
10//
11// author: Thorsten Glebe
12// HERA-B Collaboration
13// Max-Planck-Institut fuer Kernphysik
14// Saupfercheckweg 1
15// 69117 Heidelberg
16// Germany
17// E-mail: T.Glebe@mpi-hd.mpg.de
18//
19// Description: Track class
20//
21// changes:
22// 21 Aug 2000 (TG) creation
23// 24 Aug 2000 (TG) added ArtePointer declaration
24// 04 Sep 2000 (TG) added pt(), theta(), phi(), eta() member function
25// 05 Sep 2000 (TG) declared some member functions inline
26// 06 Sep 2000 (TG) added xvec(), tvec(), pvec(), E() member functions
27// 06 Sep 2000 (TG) added comments
28// 12 Sep 2000 (TG) added ex(), ey(), ez(), evec(), rm() member functions
29// 10 Okt 2000 (TG) data member t_GM and member function GM() added
30// 26 Okt 2000 (TG) added #include <cmath>
31// 26 Okt 2000 (TG) #ifndef NO_ARTE added (for integration into CLUE)
32// 14 Nov 2000 (TG) added constructors for CLUE (RecoSegment, ArtePointer<RecoTrack>)
33// 25 Nov 2000 (TG) added RecoTrack constructor
34// 24 Feb 2001 (TG) operator= added, added restmass as default CTOR argument
35// 04 Sep 2001 (TG) added TrackIf CTOR
36// 17 Okt 2001 (TG) added/improved #defines for ARTE
37// 05 Okt 2001 (TG) bug fix: added CINTOBJECT to get cintdict compilation right
38// 09 Jan 2002 (TG) changed return type of propagate() to bool
39//
40// *****************************************************************************
41
42#include <iosfwd>
43#include "vt++/CMatrix.hh"
44#include "vt++/VtVector.hh"
45#include "vt++/VtRelationList.hh"
46#include "TObject.h"
47
48#if defined USE_ROOT
49#ifndef __CINT__
50#include "smatrix/SMatrix.hh"
51#include "smatrix/SVector.hh"
52#endif
53#endif
54
55namespace VERTEX {
56 class Kalman;
57
60 //============================================================================
61 // Class Track
62 //============================================================================
63 class Track: public RelationList, public TObject
64 {
65 public:
68 Track();
70 Track(const Track* const track);
72 Track(const Track& rhs);
76 virtual ~Track();
78 void set( double x,double y,double z,double tx,double ty,double p,
79 const MATRIX::CMatrix& c);
80 Track& operator=(const Track& rhs);
81
84 float x() const;
86 float y() const;
88 float x(float z) const;
90 float y(float z) const;
92 float z() const;
94 float tx() const;
96 float ty() const;
98 float p() const;
100 float pt() const;
102 float pz() const;
104 float chi2() const;
106 unsigned short int ndf() const;
108 float phi() const;
110 float theta() const;
112 float eta() const;
114 int charge() const;
116 float energy(double rm = 0.) const;
118 float xf(double rm = 0.) const;
120 float rap(double rm = 0.) const;
121
122#if defined USE_ROOT && !defined CINTOBJECT
124 SVector<double,3> xvec() const;
126 SVector<double,3> tvec() const;
128 SVector<double,3> evec() const;
130 SVector<double,3> pvec() const;
132 const SMatrix<double,5>& COV() const;
134 const SMatrix<double,5>& CINV() const;
135#else
137 inline MATRIX::VtVector evec() const { return MATRIX::VtVector(ex(),ey(),ez()); }
139 inline MATRIX::VtVector tvec() const { return MATRIX::VtVector(tx(),ty(),1.); }
141 inline MATRIX::VtVector pvec() const { return MATRIX::VtVector(px(),py(),pz()); }
142#endif
143
145 float cov_x(double dz=0.) const;
147 float cov_y(double dz=0.) const;
149 float cov_tx() const;
151 float cov_ty() const;
153 float cov_p() const;
154
156 bool isValid() const;
158 void valid();
160 void invalid();
161
162#if defined USE_ROOT
164 bool operator==(const Track& rhs) const;
165
167// void collect(vector<Track*>& c) const;
168#endif
170 bool propagate(const double zz);
171
174 double px() const;
176 double py() const;
178 double ex() const;
180 double ey() const;
182 double ez() const;
183
185 double xerr(double dz = 0) const;
187 double yerr(double dz = 0) const;
189 double txerr() const;
191 double tyerr() const;
193 double perr() const;
194
195
197 const MATRIX::CMatrix& V() const;
199 const MATRIX::CMatrix& G() const;
201 const MATRIX::CMatrix& GM() const;
202
205 std::ostream& print( std::ostream& os ) const;
207 double rm() const;
209 void rm(const double mass);
211 double rmCC() const;
213 void rmCC(const double mass);
214
215 void delete_mom();
216
219 //bool operator< (const Track& rhs);
221 //int operator==(const Track& rhs);
222
223 private:
224 short int t_Q; // charge
225 double t_rm; // rest-mass (in Vt: xmass)
226 double t_rmCC; // a conjugated mass hypothesis
227 protected:
228 MATRIX::VtVector t_p; // track parameter vector (x,y,z,tx,ty,p)
229 // (in Vt: p(1,...,5)
230 MATRIX::CMatrix t_V; // 5x5 covariance matrix of Track
231 MATRIX::CMatrix t_G; // inverse of t_cov (5x5)
232 MATRIX::CMatrix t_GM; // inverse of t_cov (4x4, without momentum)
233
235
236 }; // class Track
237
238
239 //==============================================================================
240 // operator <<
241 //==============================================================================
242 inline std::ostream& operator<< ( std::ostream& os, const Track& t ) {
243 return t.print(os);
244 }
245
246} // end of namespace VERTEX
247
248#ifndef __CINT__
249#include "VtTrack.icc"
250#endif
251#endif
brick dz
Definition: RecDispMC.C:107
Definition: CMatrix.hh:63
Definition: VtVector.hh:45
Definition: VtRelationList.hh:227
bool operator==(const RelationList &rhs) const
Definition: VtRelationList.hh:283
Definition: VtTrack.hh:64
float eta() const
rapidity $-\log\tan(\theta/2.)$
Definition: VtTrack.C:162
float x() const
x at z=0 (in Vt: p(3,..))
Definition: VtTrack.C:148
double t_rmCC
Definition: VtTrack.hh:226
unsigned short int ndf() const
dummy function: always return 0
Definition: VtTrack.C:159
double ez() const
$1/\sqrt{1+t_x^2+t_y^2}$
Track()
Definition: VtTrack.C:64
float z() const
z = 0 in Vt
Definition: VtTrack.C:150
double perr() const
$\sqrt{\sigma^2_p}$
bool isValid() const
dummy function: always return true
Definition: VtTrack.C:195
Track(const MATRIX::VtVector &v, const MATRIX::CMatrix &c)
float cov_x(double dz=0.) const
get $\sigma_x^2|_{z+dz}$
Definition: VtTrack.C:184
const MATRIX::CMatrix & GM() const
inverse $4\times4$ covariance matrix, without momentum
Definition: VtTrack.C:207
double t_rm
Definition: VtTrack.hh:225
MATRIX::VtVector t_p
Definition: VtTrack.hh:228
MATRIX::VtVector tvec() const
$\vec{v} = (t_x,t_y,1.)$ refitted slope vector
Definition: VtTrack.hh:139
double ey() const
$t_y \cdot e_z$
MATRIX::VtVector evec() const
$\vec{v} = (e_x,e_y,e_z)$ unit vector along refitted track
Definition: VtTrack.hh:137
double ex() const
$t_x \cdot e_z$
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
virtual ~Track()
Definition: VtTrack.C:95
std::ostream & print(std::ostream &os) const
called by cout
Definition: VtTrack.C:213
void delete_mom()
Definition: VtTrack.C:244
ClassDef(Track, 1)
float cov_y(double dz=0.) const
get $\sigma_y^2|_{z+dz}$
Definition: VtTrack.C:188
double py() const
$p \cdot e_y$
void rm(const double mass)
set rest-mass (needed for mass constrained fits)
void rmCC(const double mass)
set conjugated rest-mass
MATRIX::CMatrix t_V
Definition: VtTrack.hh:230
double tyerr() const
$\sqrt{\sigma^2_{ty}}$
float phi() const
azimuthal angle $\phi$
Definition: VtTrack.C:160
double yerr(double dz=0) const
$\sqrt{\sigma_y^2(z+dz)}$
double xerr(double dz=0) const
$\sqrt{\sigma_x^2(z+dz)}$
void set(double x, double y, double z, double tx, double ty, double p, const MATRIX::CMatrix &c)
Definition: VtTrack.C:111
int charge() const
charge
Definition: VtTrack.C:163
bool propagate(const double zz)
propagate track to $z$
Definition: VtTrack.C:230
double rmCC() const
conjugated rest-mass
float energy(double rm=0.) const
Energy with given rest-mass (in GeV) $E = \sqrt{p^2 + m^2}$.
Definition: VtTrack.C:164
void invalid()
dummy function: do nothing
Definition: VtTrack.C:197
double rm() const
rest-mass
double px() const
$p \cdot e_x$
void valid()
dummy function: do nothing
Definition: VtTrack.C:196
float rap(double rm=0.) const
Rapidity $y = \frac{1}{2}\ln(\frac{E+p_z}{E-p_z})$.
Definition: VtTrack.C:166
double txerr() const
$\sqrt{\sigma^2_{tx}}$
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 theta() const
polar angle $\theta = \cos^{-1}(e_z)$
Definition: VtTrack.C:161
MATRIX::CMatrix t_G
Definition: VtTrack.hh:231
float chi2() const
dummy function: always return 0
Definition: VtTrack.C:158
float cov_ty() const
get $\sigma_{t_y}^2$
Definition: VtTrack.C:193
MATRIX::CMatrix t_GM
Definition: VtTrack.hh:232
MATRIX::VtVector pvec() const
$\vec{v} = (p_x,p_y,p_z)$ refitted mom. vector
Definition: VtTrack.hh:141
float pz() const
$p \cdot e_z$
Definition: VtTrack.C:157
Track & operator=(const Track &rhs)
Definition: VtTrack.C:128
short int t_Q
Definition: VtTrack.hh:224
float cov_p() const
get $\sigma_p^2$
Definition: VtTrack.C:194
float pt() const
transv. momentum $\sqrt{p_x^2 + p_y^2}$
Definition: VtTrack.C:156
float xf(double rm=0.) const
Feynman Variable $x_F = \frac{E+p_z}{(E+p_z)_{max}}$.
Definition: VtTrack.C:165
float tx() const
slope (in Vt: p(1,..))
Definition: VtTrack.C:153
float cov_tx() const
get $\sigma_{t_x}^2$
Definition: VtTrack.C:192
Definition: bitview.h:14
TTree * t
Definition: check_shower.C:4
float mass
Definition: check_vertex.C:21
Definition: VtDistance.hh:30
std::ostream & operator<<(std::ostream &os, const VtIni &t)
Definition: VtIni.hh:83