FEDRA emulsion software from the OPERA Collaboration
SVertex.hh
Go to the documentation of this file.
1#ifndef __SVERTEX_HH
2#define __SVERTEX_HH
3// ********************************************************************
4//
5// source:
6//
7// type: source code
8//
9// created: 06. Apr 2001
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: A fixed size Vertex class
20//
21// changes:
22// 06 Apr 2001 (TG) creation
23// 07 Apr 2001 (TG) mass() added
24// 10 Apr 2001 (TG) findVertexVt(), kalman_, v_bk13, template parameter use_mom
25// added
26// 12 Apr 2001 (TG) set_track, track(), calc_mother(), calc_mother_tr(),
27// bigcov(), mother() added
28// 17 Apr 2001 (TG) Track inheritance added, calc_mother_cov() added
29// 18 Apr 2001 (TG) removed use_mom template parameter, added mass_tr
30// 24 Apr 2001 (TG) valid_, valid(), invalid() added
31// 24 Apr 2001 (TG) renamed calc_mother() to calc_mother_tr() and
32// calc_mother_tr() to calc_mother_trtr(). Added new calc_mother().
33// 24 Apr 2001 (TG) added prob()
34// 30 Apr 2001 (TG) added pt() and xfabs()
35// 02 Mai 2001 (TG) phi(), theta(), eta() added
36// 03 Mai 2001 (TG) added read/write vposR(), added findVertexNe()
37// 03 Mai 2001 (TG) CTOR added, size added
38// 08 Jun 2001 (TG) added massError()
39// 11 Jun 2001 (TG) added cov_tx(), cov_ty(), cov_p()
40// 28 Jun 2001 (TG) added charge(), pz(), energy(), xf(), rap(), removed xfabs()
41// 04 Jul 2001 (TG) moth_tr, moth_cov added
42// 05 Jul 2001 (TG) collect() added
43// 09 Jul 2001 (TG) added valid(), invalid(), isValid() from Track, removed
44// valid(), invalid() from Vertex, added ntracks()
45// 10 Jul 2001 (TG) operator==() added
46// 23 Aug 2001 (TG) added CTOR, added copy CTOR
47// 03 Sep 2001 (TG) added operator=()
48// 02 Okt 2001 (TG) set_vpos() added
49// 07 Nov 2001 (TG) documentation improved
50// 03 Jan 2002 (TG) isMotherTr(), isMotherCov() added, renamed valid_ to
51// validKalmanFit, changed return type of calcVertex2D
52// 04 Jan 2002 (TG) added isKalmanFit()
53// 09 Jan 2002 (TG) changed return type of propagate() to bool
54// 11 Jan 2002 (TG) added v_ndf member
55// 16 Jan 2002 (TG) added operator=(SVertex<NTR>)
56//
57// ********************************************************************
58#include <cmath>
59#include <iosfwd>
60#include "vt++/VtVertex.hh"
61#include "vt++/VtTrack.hh"
62#include "smatrix/SVector.hh"
63#include "smatrix/SKalman.hh"
64
69//==============================================================================
70// SVertex
71//==============================================================================
72template <unsigned int NTR>
73class SVertex : public Vertex, public Track {
74public:
79 SVertex(const SVertex<NTR>& rhs);
81 SVertex(const Track& t);
83 SVertex(const Track& t1, const Track& t2);
85 SVertex(const Track& t1, const Track& t2, const Track& t3);
87 SVertex(const Track& t1, const Track& t2, const Track& t3,
88 const Track& t4);
90 SVertex(const SVertex<NTR-1>& vtx, const Track& t);
91
94
97 static const unsigned int size = NTR;
98
101 float vx() const;
103 float vy() const;
105 float vz() const;
107 float chi2() const;
109 unsigned short int ndf() const;
111 float prob() const;
113 unsigned short int ntracks() const;
115 const Track* track(unsigned int i) const;
117 const Track*& track(unsigned int i);
118
122 const SMatrix<double,3>& VCOV() const;
123
125 float vtx_cov_x() const;
127 float vtx_cov_y() const;
129 float vtx_cov_z() const;
130
133 const SVector<double,3>& vposR() const;
137 void set_vpos(const SVector<double,3>& pos);
138
140 void set_track(unsigned int i, const Track& t);
142 const SKalman<NTR>& kalman(unsigned int i) const;
144 SKalman<NTR>& kalman(unsigned int i);
146 const SMatrix<double,3>& VCINV() const;
148 const SVector<double,6>& mother() const;
151
168 double mass(const SVector<double,NTR>& rm) const;
170 double massError(const SVector<double,NTR>& rm) const;
172 double mass_tr(const SVector<double,NTR>& rm) const;
173
179 float x() const;
181 float y() const;
183 float x(float z) const;
185 float y(float z) const;
187 float z() const;
189 float tx() const;
191 float ty() const;
193 float p() const;
195 float pt() const;
197 float pz() const;
199 float phi() const;
201 float theta() const;
203 float eta() const;
205 int charge() const;
206
208 float energy(double mass = 0.) const;
211 float xf(double mass = 0.) const;
213 float rap(double mass = 0.) const;
214
216 float cov_x(double dz=0.) const;
218 float cov_y(double dz=0.) const;
220 float cov_tx() const;
222 float cov_ty() const;
224 float cov_p() const;
225
230 bool isValid() const;
232 void valid();
234 void invalid();
235
244
246 const SMatrix<double,5>& COV() const;
248 const SMatrix<double,5>& CINV() const;
249
251 void collect(vector<Track*>& c) const;
253 bool propagate(const double z);
255 std::ostream& print( std::ostream& ) const;
256
259 bool operator==(const Track& rhs) const;
260
263 bool isMotherTr() const;
265 bool isMotherCov() const;
267 bool isKalmanFit() const;
268
270 bool filter();
272 bool smooth();
274 bool smoothC();
276 inline double bk13() const { return v_bk13; }
286 bool bigcov();
287
288private:
289 SVector<double,3> vpos_; // vertex position
290 double v_bk13; // vtx-z before kalman filter could be vz?
291 const Track* tracks[NTR]; // track pointers
292 SKalman<NTR> kalman_[NTR]; // kalman objects
293 SMatrix<double,3> v_CS; // vertex cov. matrix
294 SMatrix<double,3> v_CINV; // inverse vertex cov. matrix
295 double v_chi2; // chi2 of vertex fit
296 unsigned int v_ndf; // no of degrees of vertex fit
297 SVector<double,6> t_p; // track parameter of mother particle
298 SMatrix<double,3*NTR+3> v_covn; // cov. matrix of all track parameters
299 SMatrix<double,5> cov_; // cov. matrix of reconstructed track
300 SMatrix<double,5> cinv_; // inv. cov. matrix of reconstructed track
301 bool validKalmanFit; // was Kalman fit successful?
302 bool valid_; // is vertex valid?
303 bool moth_tr; // is mother track calculated?
304 bool moth_cov; // is mother covariance matrix calculated?
305 bool validity; // vertex validity, manually set
306}; // class SVertex
307
308
309//==============================================================================
310// operator <<
311//==============================================================================
312template <unsigned int NTR>
313inline std::ostream& operator<< (std::ostream& os, const SVertex<NTR>& vt) {
314 return vt.print( os );
315}
316
317#include "SVertex.icc"
318#endif
std::ostream & operator<<(std::ostream &os, const SVertex< NTR > &vt)
Definition: SVertex.hh:313
Definition: SKalman.hh:49
Definition: SVertex.hh:73
SVector< double, 3 > evec() const
$\vec{e} = (e_x,e_y,ez)$ unit vector along $\vec{p}$ of mother track
float theta() const
polar angle $\theta = \cos^{-1}(e_z)$ [deg]
SVector< double, 3 > vpos_
Definition: SVertex.hh:289
bool valid_
Definition: SVertex.hh:302
const SMatrix< double, 3 > & VCOV() const
vertex cov. matrix
float pz() const
mother track momentum along
bool isValid() const
float cov_x(double dz=0.) const
get $\sigma_x^2|_{z+dz}$
SVector< double, 6 > t_p
Definition: SVertex.hh:297
float ty() const
mother track slope $t_y$
SMatrix< double, 3 *NTR+3 > v_covn
Definition: SVertex.hh:298
float y() const
$y$ of reconstructed track
float p() const
mother track momentum $p$
float chi2() const
vertex $\chi^2$
bool filter()
do a filter step
float cov_ty() const
get $\sigma_{t_y}^2$
float eta() const
rapidity $\eta = -\log\tan(\theta/2.)$
const SKalman< NTR > & kalman(unsigned int i) const
read only access to Kalman objects
float prob() const
upper tail $\chi^2$ probability
unsigned short int ndf() const
degrees of freedom of vertex fit
SVertex< NTR > & operator=(const SVertex< NTR > &rhs)
bool findVertex3D(SVector< double, 3 > &pos) const
return vertex computed with 3D analytical method
SMatrix< double, 5 > cinv_
Definition: SVertex.hh:300
int charge() const
track charge (=sum of track charges in vertex)
float z() const
$z$ of reconstructed track
unsigned short int ntracks() const
no of tracks in vertex
float vtx_cov_z() const
$\sigma_z^2$ of vertex
const Track *& track(unsigned int i)
read/write track access
SVector< double, 3 > findVertex2D() const
return vertex computed with 2D analytical method
const Track * tracks[NTR]
Definition: SVertex.hh:291
bool operator==(const Track &rhs) const
compare Track pointers
const SMatrix< double, 3 *NTR+3 > & covn() const
for internal use only
void collect(vector< Track * > &c) const
collect pointers
const SMatrix< double, 3 > & VCINV() const
inverse vertex cov. matrix
SMatrix< double, 3 > v_CS
Definition: SVertex.hh:293
const SMatrix< double, 5 > & COV() const
covariance matrix of mother track
float vx() const
vertex $x$ position $v_x$
SVertex(const Track &t)
create Vertex object & set mother track by $t$
const SVector< double, 6 > & mother() const
mother track parameters $\vec{m}=(x,y,z,t_x,t_y,p)$
const SVector< double, 3 > & vposR() const
vertex position $\vec{v} = (v_x,v_y,v_z)$ (fast readonly access)
unsigned int v_ndf
Definition: SVertex.hh:296
bool moth_cov
Definition: SVertex.hh:304
void set_vpos(const SVector< double, 3 > &pos)
set vertex position manually
void valid()
mark mother track as valid
SVector< double, 3 > EstimateVertex() const
compute vertex with error weighted 2D analytical method
float tx() const
mother track rack slope $t_x$
double mass_tr(const SVector< double, NTR > &rm) const
compute invariant mass using measured track momenta
SVector< double, 3 > vpos() const
vertex position $\vec{v} = (v_x,v_y,v_z)$
SVector< double, 3 > pvec() const
mother track momentum vector $\vec{p} = (p_x,p_y,p_z)$
void set_track(unsigned int i, const Track &t)
set track list
double v_chi2
Definition: SVertex.hh:295
bool findVertexNe()
Vt based Kalman filter vertex fit, without first estimation.
bool isMotherCov() const
is mother cov. matrix calculated?
std::ostream & print(std::ostream &) const
used by operator<<()
bool calc_mother()
calculate mother track + cov. matrix
float rap(double mass=0.) const
Rapidity $y = \frac{1}{2}\ln(\frac{E+p_z}{E-p_z})$.
void invalid()
mark mother track as invalid
SVector< double, 3 > xvec() const
position $\vec{x} = (x,y,z)$ of mother track
SVertex(const SVertex< NTR-1 > &vtx, const Track &t)
add a track
SVertex(const SVertex< NTR > &rhs)
float cov_y(double dz=0.) const
get $\sigma_y^2|_{z+dz}$
float cov_p() const
get $\sigma_p^2$
const Track * track(unsigned int i) const
read only track access
float vy() const
vertex $y$ position $v_y$
bool validKalmanFit
Definition: SVertex.hh:301
bool smoothC()
needed for bigcov()
SKalman< NTR > kalman_[NTR]
Definition: SVertex.hh:292
bool calcVertex2D()
compute vertex with a 2D analytical method
double v_bk13
Definition: SVertex.hh:290
SVertex(const Track &t1, const Track &t2, const Track &t3, const Track &t4)
bool findVertexVt()
Vt based Kalman filter vertex fit.
float vtx_cov_y() const
$\sigma_y^2$ of vertex
SVertex(const Track &t1, const Track &t2, const Track &t3)
float cov_tx() const
get $\sigma_{t_x}^2$
float x(float z) const
linear extrapolation of $x$ to z
static const unsigned int size
Definition: SVertex.hh:97
double bk13() const
Definition: SVertex.hh:276
float y(float z) const
linear extrapolation of $y$ to z
bool propagate(const double z)
propagate reconsructed track to $z$
SVector< double, 3 > tvec() const
slope $\vec{t} = (t_x,t_y,1)$ of mother track
bool isKalmanFit() const
has Kalman fit succeeded?
bool calc_mother_trtr()
calculate mother track using track info
SVector< double, 2 > rmsDistAngle() const
return rms distance and rms opening angle of tracks in vertex
float vz() const
vertex $z$ position $v_z$
bool validity
Definition: SVertex.hh:305
float x() const
$x$ of reconstructed track
bool isMotherTr() const
is mother track calculated?
SMatrix< double, 5 > cov_
Definition: SVertex.hh:299
SKalman< NTR > & kalman(unsigned int i)
read/write access to Kalman objects
float pt() const
transverse mother track momentum $p_t$
float energy(double mass=0.) const
$E = \sqrt{m^2 + p^2}$
const SMatrix< double, 5 > & CINV() const
inverse cov. matrix of mother track
bool calc_mother_tr()
calculate mother track using kalman info
float xf(double mass=0.) const
float phi() const
azimuthal angle $\phi$ [deg]
float vtx_cov_x() const
$\sigma_x^2$ of vertex
bool moth_tr
Definition: SVertex.hh:303
bool calc_mother_cov()
calculate covariance matrix of mother track
double massError(const SVector< double, NTR > &rm) const
compute invariant mass error caused by momentum error
double mass(const SVector< double, NTR > &rm) const
compute invariant mass using refitted track momenta
SVertex(const Track &t1, const Track &t2)
SVector< double, 3 > & vposR()
vertex position $\vec{v} = (v_x,v_y,v_z)$ (fast read/write access)
bool bigcov()
construct $(3\cdot n+3)\times(3\cdot n+3)$ cov. matrix
bool smooth()
smooth vertex parameters
SMatrix< double, 3 > v_CINV
Definition: SVertex.hh:294
Definition: Track.h:10
double dz
Definition: Track.h:46
TTree * t
Definition: check_shower.C:4