FEDRA emulsion software from the OPERA Collaboration
VtVertex.hh
Go to the documentation of this file.
1#ifndef __VTVERTEX_HH
2#define __VTVERTEX_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: Vertex class
20//
21// changes:
22// 21 Aug 2000 (TG) creation
23// 25 Aug 2000 (TG) remove_worst added
24// 04 Sep 2000 (TG) member function ndf() added. chi2() returns now chi2
25// and not chi2/ndf any more
26// 04 Sep 2000 (TG) Modified argument of VtRemoveTrack()
27// 05 Sep 2000 (TG) declared some member functions as inline
28// 05 Sep 2000 (TG) erase() member function added
29// 12 Sep 2000 (TG) private member function bigdim() added
30// 12 Sep 2000 (TG) member function bigcov() added
31// 13 Sep 2000 (TG) added xmass member and add() member function
32// 14 Sep 2000 (TG) added v_CINV member
33// 14 Sep 2000 (TG) bug fix in documentation
34// 02 Okt 2000 (TG) added chi2l() member function
35// 04 Okt 2000 (TG) removed v_use_momentum data member
36// 17 Okt 2000 (TG) modified comments
37// 18 Okt 2000 (TG) MassConstr() member function added
38// 24 Okt 2000 (TG) invalid() member function added, v_covn, v_Mvalid member added,
39// changed return type of bigcov()
40// 25 Okt 2000 (TG) covn() member function added
41// 26 Okt 2000 (TG) #ifndef NO_ARTE added (for integration into CLUE)
42// 14 Nov 2000 (TG) added constructors for CLUE (RecoSegment, ArtePointer<RecoTrack>)
43// 31 Jan 2001 (TG) added prob() member function
44// 02 Feb 2001 (TG) added global ClMember declarations
45// 20 Feb 2001 (TG) pmaxfrac() added, ClMember pmaxfrac added
46// 27 Feb 2001 (TG) improved ndf() (mass constrained vertex fit case added)
47// 28 Feb 2001 (TG) added default argument to calc_mother()
48// 06 Mar 2001 (TG) v_angdist added, angle() and dist() check now whether they
49// are computed.
50// 12 Mar 2001 (TG) set_invalid() added
51// 09 Jul 2001 (TG) changed return type of ntracks() to unsigned int, added clear()
52// 04 Sep 2001 (TG) added push_back(TrackIf), added VertexIf inheritance
53// 17 Okt 2001 (TG) added/improved #defines for ARTE
54//
55// *****************************************************************************
56
57#include <iosfwd>
58#include <vector>
59#include "CMatrix.hh"
60#include "VtVector.hh"
61#include "VtTrack.hh"
62#include "VtMassC_list.hh"
63
64#ifdef WIN32
65# include "Rtypes.h"
66#endif
67
68void Dinv1(int *n, double *a, int *idim, double *ir, int *ifail);
69
70namespace VERTEX {
71 //
72 // Some Declarations
73 //
74 class Track;
75 class Vertex;
76
77 typedef std::vector<double> Vector_d;
78
79 typedef std::vector<Track*> Track_v;
82
84 //============================================================================
85 // Class Vertex
86 //============================================================================
87 class Vertex: public Track
88 {
89 public:
90
93 Vertex();
95 Vertex(const Vertex& rhs);
97 Vertex(const Track_v& v);
99 Vertex(Track& t1, Track& t2);
101 ~Vertex();
102
105 float vx() const;
107 float vy() const;
109 float vz() const;
111 float chi2() const;
113 float prob() const;
115 unsigned short int ndf() const;
117 unsigned short int ntracks() const;
118
119#if defined USE_ROOT && !defined __CINT__
121 SVector<double,3> vpos() const;
123 const SMatrix<double,3>& VCOV() const;
124#endif
125
127 float vtx_cov_x() const;
129 float vtx_cov_y() const;
131 float vtx_cov_z() const;
132
136 double angle() const;
138 double dist() const;
140 double vxerr() const;
142 double vyerr() const;
144 double vzerr() const;
145
147 double pmaxfrac() const;
148
150 const MATRIX::VtSymMatrix& CS() const;
152 const MATRIX::VtSymMatrix& covn() const;
154 bool valid() const;
157
162 double mass(const bool use = false) const;
164 double massCC(const bool use = false) const;
166 double mass(double m1) const;
168 double mass(double m1, double m2) const;
170 double mass(double m1, double m2, double m3) const;
172 double mass(double m1, double m2, double m3, double m4) const;
174 double mass(const Vector_d& m02, bool use = false, bool CC = false) const;
175
178 const bool findVertex2D(void);
180 const bool findVertex3D(void);
182 const bool findVertexVt ();
184 const bool VtEstimateVertex();
186 const bool VtEstimateVertexMath(double& x, double& y, double& z);
188 const bool VtEstimateVertexMathTA(double& x, double& y, double& z);
190 const bool VtMass();
191
194 void push_back(Track& track);
195
199 double rmsDistAngle() const;
201 bool calc_mother(bool use_kalman = true);
203 bool calc_mother_cov();
204
206 const bool remove_last();
208 void remove_worst();
210 double track_chi2(int i);
212 int track_worst();
214 void use_momentum(const bool use);
216 void use_kalman(const bool use);
218 bool use_kalman() const;
220 double distance(double x, double y, double z) const;
222 double distance(Track& t) const;
224 double distance(const Vertex& rhs) const;
225
228 void clear();
230 const bool VtFilter ();
232 const bool VtInverseFilter () const;
234 void VtSmoothX ();
236 void VtSmoothQ ();
238 const bool VtRemoveTrack (Relation& track);
240 const bool VtFit ();
242 const MATRIX::VtVector& xv() const;
244 const MATRIX::VtVector& xvs() const;
246 std::ostream& print(std::ostream& os) const;
248 double chi2n() const;
250 double chi2l() const;
252 virtual void remove(Relation* r);
254 virtual const iterator erase(const iterator& pos);
257
260 void MassConstr(double m);
262 MassC& addMassConstr(double m=0.);
264 void clearMassConstr();
266 unsigned int nMassConstr() const;
267
270 const Vertex& operator= (const Vertex& rhs);
272 const bool operator==(const Vertex& rhs) const;
274 const double operator- (const Vertex& rhs) const;
275
276 void add_track(Track &t);
277 const Track *get_track(int i) const;
278
279 private:
280 const unsigned int bigdim() const;
281 // flag vertex as invalid
282 void invalid();
283
284 bool v_use_kalman; // use refitted track parameters
285 bool v_valid; // is vertex information valid?
286 bool v_Mvalid; // is mass constrained vertex information valid?
287 bool v_mother; // is the mother particle calculated?
288 mutable bool v_angdist; // is angle and dist calculated?
289 mutable double v_angle; // rms of opening angle
290 mutable double v_dist; // rms of distance to vertex
291 double v_bk13; // vtx-z at begin of iteration step
292 // (in Vt: Bk13(iv), b13)
293 double v_chi2; // chi2 of vertex-fit
294 MATRIX::VtVector kal_xv; // x,y,z of vtx from kalman filter
295 MATRIX::VtVector kal_xvs; // pos of Vertex after smooth
296 // (xv-smoothed) (in Vt: xvs(1..3,...))
297 MATRIX::VtSymMatrix v_CS; // v_C smoothed (=v_C of last
298 // filtered track) (in Vt: CS)
299 MATRIX::VtSymMatrix v_CINV; // inverse of CS at initialization time
300 MATRIX::VtSymMatrix* v_covn; // (3*n+3)x(3*n+3) big cov. matrix
301 Track_v tracks; // tracks owned by vertex
302 MassC_v xmass; // vector of mass constraint pointer
303
304#ifdef WIN32
305 ClassDef(Vertex, 0); //GS comment if <VersionNumber> >= 1 crash with Windows Vista
306#else
308#endif
309 }; // class Vertex
310
311
312 //==============================================================================
313 // operator<<
314 //==============================================================================
315 inline std::ostream& operator<<(std::ostream& s, const Vertex& a) {
316 return a.print(s);
317 }
318
319} // namespace VERTEX
320
321#ifndef __CINT__
322#include "VtVertex.icc"
323#endif
324
325#endif
void Dinv1(int *n, double *a, int *idim, double *ir, int *ifail)
Definition: VtSqMatrix.C:321
void a()
Definition: check_aligned.C:59
Definition: VtSymMatrix.hh:49
Definition: VtVector.hh:45
Definition: Track.h:10
Definition: VtMassC.hh:50
Definition: VtRelationList.hh:67
Definition: VtRelation.hh:51
Definition: VtTrack.hh:64
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
float y() const
y at z=0 (in Vt: p(4,..))
Definition: VtTrack.C:149
Definition: VtVertex.hh:88
void use_kalman(const bool use)
use refitted track parameters or not
double pmaxfrac() const
$\max|p_i/p_j|$ max. momentum fraction
Definition: VtVertex.C:1570
const bool remove_last()
remove last track from vertex, do inverse Kalman filter
Definition: VtVertex.C:725
const bool VtEstimateVertex()
weighted 2-D min. dist.
Definition: VtVertex.C:1241
double massCC(const bool use=false) const
rest-mass 0 or conjugated track rest mass
Definition: VtVertex.C:925
const MATRIX::VtSymMatrix & bigcov()
compute general $(3\cdot n+3) \times (3\cdot n+3)$ covariance matrix
Definition: VtVertex.C:341
double vyerr() const
$\sqrt{\sigma_{vy}^2}$ vertex $y$-error
bool v_mother
Definition: VtVertex.hh:287
const bool operator==(const Vertex &rhs) const
Definition: VtVertex.C:1555
const MATRIX::VtVector & xv() const
MATRIX::VtSymMatrix * v_covn
Definition: VtVertex.hh:300
double track_chi2(int i)
get track $\chi^2$ contribution
Definition: VtVertex.C:783
void clearMassConstr()
delete the mass constraints
Definition: VtVertex.C:328
unsigned short int ndf() const
degrees of freedom of vertex fit
Definition: VtVertex.C:238
unsigned int nMassConstr() const
number of mass constraints
float chi2() const
$\chi^2$ of vertex fit
Definition: VtVertex.C:236
const bool findVertex3D(void)
unweighted 3-D min. dist.
Definition: VtVertex.C:1011
void add_track(Track &t)
Definition: VtVertex.C:261
const bool VtEstimateVertexMath(double &x, double &y, double &z)
estimate Vertex without changing Vertex object
Definition: VtVertex.C:1305
MATRIX::VtVector kal_xvs
Definition: VtVertex.hh:295
const bool VtInverseFilter() const
Definition: VtVertex.C:1419
const double operator-(const Vertex &rhs) const
Definition: VtVertex.C:1563
bool v_angdist
Definition: VtVertex.hh:288
bool v_Mvalid
Definition: VtVertex.hh:286
double v_angle
Definition: VtVertex.hh:289
bool v_valid
Definition: VtVertex.hh:285
MassC_v xmass
Definition: VtVertex.hh:302
double distance(double x, double y, double z) const
$\chi^2$ distance to space point, $ndf = 3$
Definition: VtVertex.C:803
int track_worst()
get track with biggest $\chi^2$ contribution
Definition: VtVertex.C:758
double vzerr() const
$\sqrt{\sigma_{vz}^2}$ vertex $z$-error
const MATRIX::VtVector & xvs() const
vector of vertex-positions
MassC & addMassConstr(double m=0.)
add a mass constraint
Definition: VtVertex.C:320
float vtx_cov_z() const
$\sigma_z^2$ of vertex
Definition: VtVertex.C:256
float vz() const
$z$ of vertex
Definition: VtVertex.C:235
double angle() const
double chi2l() const
$\chi^2$ contribution of last filter step
Definition: VtVertex.C:865
float vy() const
$y$ of vertex
Definition: VtVertex.C:234
MATRIX::VtVector kal_xv
Definition: VtVertex.hh:294
const bool findVertex2D(void)
unweighted 2-D min. dist.
Definition: VtVertex.C:1111
const bool VtFit()
Kalman filter + smoother.
Definition: VtVertex.C:1162
void use_momentum(const bool use)
set for all tracks whether momentum should be used or not
Definition: VtVertex.C:279
void push_back(Track &track)
Definition: VtVertex.C:271
const bool VtFilter()
Definition: VtVertex.C:1378
Vertex(const Track_v &v)
double mass(const bool use=false) const
rest-mass 0 or track rest mass
Definition: VtVertex.C:929
unsigned short int ntracks() const
no of tracks in vertex $n$
Definition: VtVertex.C:240
float vtx_cov_x() const
$\sigma_x^2$ of vertex
Definition: VtVertex.C:254
const bool VtMass()
Kalman filter with mass constraints.
Definition: VtVertex.C:535
double chi2n() const
$\chi^2$ after last filter step
Definition: VtVertex.C:858
double v_bk13
Definition: VtVertex.hh:291
void MassConstr(double m)
add a mass constraint for all tracks, $m$ = mass of mother particle
Definition: VtVertex.C:306
float vtx_cov_y() const
$\sigma_y^2$ of vertex
Definition: VtVertex.C:255
const unsigned int bigdim() const
const Track * get_track(int i) const
Definition: VtVertex.C:266
std::ostream & print(std::ostream &os) const
called by operator<<()
Definition: VtVertex.C:875
float vx() const
$x$ of vertex
Definition: VtVertex.C:233
Vertex()
Definition: VtVertex.C:88
bool calc_mother(bool use_kalman=true)
calc track params of mother particle
Definition: VtVertex.C:392
double rmsDistAngle() const
calc rms dist and rms angle
Definition: VtVertex.C:1073
const MATRIX::VtSymMatrix & CS() const
$3\times3$ Vertex covariance matrix
double dist() const
available after call to rmsDistAngle()
const bool VtEstimateVertexMathTA(double &x, double &y, double &z)
estimate Vertex without changing Vertex object
Definition: VtVertex.C:1249
ClassDef(Vertex, 1)
void remove_worst()
remove track with biggest $\chi^2$ contribution, do inverse Kalman filter
Definition: VtVertex.C:732
bool valid() const
is vertex valid?
Track_v tracks
Definition: VtVertex.hh:301
bool calc_mother_cov()
calc cov. matrix of mother particle
Definition: VtVertex.C:458
float prob() const
upper tail $\chi^2$ probability
Definition: VtVertex.C:237
const Vertex & operator=(const Vertex &rhs)
Definition: VtVertex.C:1525
MATRIX::VtSymMatrix v_CINV
Definition: VtVertex.hh:299
void clear()
clear all track-vertex relations, makes vertex invalid
Definition: VtVertex.C:196
~Vertex()
Definition: VtVertex.C:178
double vxerr() const
$\sqrt{\sigma_{vx}^2}$ vertex $x$-error
void VtSmoothQ()
smooth momenta
Definition: VtVertex.C:1506
const MATRIX::VtSymMatrix & covn() const
$(3\cdot n+3)\times(3\cdot n+3)$ general covariance matrix
virtual void remove(Relation *r)
unbook relation
Definition: VtVertex.C:290
void set_invalid()
mark vertex as invalid
void VtSmoothX()
smooth vertex position
Definition: VtVertex.C:1474
double v_dist
Definition: VtVertex.hh:290
bool v_use_kalman
Definition: VtVertex.hh:284
const bool VtRemoveTrack(Relation &track)
remove track using inverse Kalman filter
Definition: VtVertex.C:1443
MATRIX::VtSymMatrix v_CS
Definition: VtVertex.hh:297
void invalid()
Definition: VtVertex.C:218
virtual const iterator erase(const iterator &pos)
erase a relation without refit
Definition: VtVertex.C:298
const bool findVertexVt()
Vt Kalman-filter.
Definition: VtVertex.C:1174
bool use_kalman() const
double v_chi2
Definition: VtVertex.hh:293
Definition: bitview.h:14
TTree * t
Definition: check_shower.C:4
s
Definition: check_shower.C:55
Definition: VtDistance.hh:30
std::vector< Track * >::iterator Track_it
Definition: VtVertex.hh:80
std::vector< Track * > Track_v
Definition: VtVertex.hh:79
std::ostream & operator<<(std::ostream &os, const VtIni &t)
Definition: VtIni.hh:83
ConstRelationIterator const_iterator
Definition: VtRelationList.hh:54
std::list< MassC * > MassC_v
Definition: VtMassC_list.hh:28
RelationIterator iterator
Definition: VtRelationList.hh:53
std::vector< Track * >::const_iterator Track_cit
Definition: VtVertex.hh:81
std::vector< double > Vector_d
Definition: VtVertex.hh:75
void r(int rid=2)
Definition: test.C:201