CMS 3D CMS Logo

Vertex.h
Go to the documentation of this file.
1 #ifndef VertexReco_Vertex_h
2 #define VertexReco_Vertex_h
3 
20 #include <Rtypes.h>
26 #include <Math/GenVector/PxPyPzE4D.h>
27 #include <Math/GenVector/PxPyPzM4D.h>
29 
30 namespace reco {
31 
32  class Track;
33 
34  class Vertex {
35  public:
37  typedef std::vector<TrackBaseRef >::const_iterator trackRef_iterator;
41  enum { dimension = 3, dimension4D = 4 };
51  enum { size = dimension * ( dimension + 1 ) / 2, size4D = ( dimension4D ) * ( dimension4D + 1 ) / 2 };
53  typedef unsigned int index;
57  Vertex(): chi2_( 0.0 ), ndof_( 0 ), position_(0.,0.,0.), time_(0.) { validity_ = false; for(int i = 0; i < size4D; ++ i ) covariance_[i]=0.;
58 }
60  Vertex( const Point &, const Error &);
62  Vertex( const Point &, const Error4D &, double);
64  Vertex( const Point &, const Error &, double chi2, double ndof, size_t size );
66  Vertex( const Point &, const Error4D &, double time, double chi2, double ndof, size_t size );
68  bool isValid() const {return validity_;}
72  bool isFake() const {return (chi2_==0 && ndof_==0 && tracks_.empty());}
74  void add( const TrackBaseRef & r, float w=1.0 );
76  void add( const TrackBaseRef & r, const Track & refTrack, float w=1.0 );
77  void removeTracks();
78 
80  template<typename TREF>
81  float trackWeight ( const TREF & r ) const {
82  int i=0;
83  for(auto const & t : tracks_) {
84  if ( (r.id()==t.id()) & (t.key()==r.key()) ) return weights_[i]/255.f;
85  ++i;
86  }
87  return 0;
88  }
90  trackRef_iterator tracks_begin() const;
92  trackRef_iterator tracks_end() const;
94  size_t tracksSize() const;
96  const TrackBaseRef& trackRefAt(size_t idx) const { return tracks_[idx]; }
98  double chi2() const { return chi2_; }
105  double ndof() const { return ndof_; }
107  double normalizedChi2() const { return ndof_ != 0 ? chi2_ / ndof_ : chi2_ * 1e6; }
109  const Point & position() const { return position_; }
111  double x() const { return position_.X(); }
113  double y() const { return position_.Y(); }
115  double z() const { return position_.Z(); }
117  double t() const { return time_; }
119  double xError() const { return sqrt( covariance(0, 0) ); }
121  double yError() const { return sqrt( covariance(1, 1) ); }
123  double zError() const { return sqrt( covariance(2, 2) ); }
125  double tError() const { return sqrt( covariance(3, 3) ); }
127  // Note that:
128  // double error( int i, int j ) const
129  // is OBSOLETE, use covariance(i, j)
130  double covariance( int i, int j ) const {
131  return covariance_[ idx( i, j ) ];
132  }
134  CovarianceMatrix covariance() const { Error m; fill( m ); return m; }
136  CovarianceMatrix4D covariance4D() const { Error4D m; fill( m ); return m; }
137 
139  Error error() const { Error m; fill( m ); return m; }
141  Error4D error4D() const { Error4D m; fill( m ); return m; }
142 
144  void fill( CovarianceMatrix & v ) const;
146  void fill( CovarianceMatrix4D & v ) const;
147 
149  bool hasRefittedTracks() const {return !refittedTracks_.empty();}
150 
153  TrackBaseRef originalTrack(const Track & refTrack) const;
154 
157  Track refittedTrack(const TrackBaseRef & track) const;
158 
161  Track refittedTrack(const TrackRef & track) const;
162 
164  const std::vector<Track> & refittedTracks() const { return refittedTracks_;}
165 
167  math::XYZTLorentzVectorD p4(float mass=0.13957018,float minWeight=0.5) const;
168 
170  unsigned int nTracks(float minWeight=0.5) const;
171 
172  class TrackEqual {
173  public:
174  TrackEqual( const Track & t) : track_( t ) { }
175  bool operator()( const Track & t ) const { return t.pt()==track_.pt();}
176  private:
177  const Track & track_;
178  };
179 
180  private:
182  float chi2_;
184  float ndof_;
186  Point position_;
190  std::vector<TrackBaseRef> tracks_;
192  std::vector<Track> refittedTracks_;
193  std::vector<uint8_t> weights_;
195  bool validity_;
196  double time_;
197 
199  index idx( index i, index j ) const {
200  int a = ( i <= j ? i : j ), b = ( i <= j ? j : i );
201  return b * ( b + 1 ) / 2 + a;
202  }
203  };
204 
205 }
206 
207 #endif
void fill(CovarianceMatrix &v) const
fill SMatrix
float chi2_
chi-sqared
Definition: Vertex.h:182
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< double > > XYZTLorentzVectorD
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:14
const double w
Definition: UKUtility.cc:23
trackRef_iterator tracks_end() const
last iterator over tracks
Definition: Vertex.cc:81
double zError() const
error on z
Definition: Vertex.h:123
TrackBaseRef originalTrack(const Track &refTrack) const
Definition: Vertex.cc:108
Track refittedTrack(const TrackBaseRef &track) const
double y() const
y coordinate
Definition: Vertex.h:113
bool isValid() const
Tells whether the vertex is valid.
Definition: Vertex.h:68
std::vector< Track > refittedTracks_
The vector of refitted tracks.
Definition: Vertex.h:192
bool hasRefittedTracks() const
Checks whether refitted tracks are stored.
Definition: Vertex.h:149
const TrackBaseRef & trackRefAt(size_t idx) const
python friendly track getting
Definition: Vertex.h:96
math::Error< dimension >::type Error
covariance error matrix (3x3)
Definition: Vertex.h:43
std::vector< TrackBaseRef > tracks_
reference to tracks
Definition: Vertex.h:190
ErrorD< N >::type type
Definition: Error.h:33
CovarianceMatrix4D covariance4D() const
return SMatrix 4D
Definition: Vertex.h:136
unsigned int index
index type
Definition: Vertex.h:53
double covariance(int i, int j) const
(i, j)-th element of error matrix, i, j = 0, ... 2
Definition: Vertex.h:130
const std::vector< Track > & refittedTracks() const
Returns the container of refitted tracks.
Definition: Vertex.h:164
const Point & position() const
position
Definition: Vertex.h:109
math::Error< dimension >::type CovarianceMatrix
covariance error matrix (3x3)
Definition: Vertex.h:45
TrackEqual(const Track &t)
Definition: Vertex.h:174
const Track & track_
Definition: Vertex.h:177
T sqrt(T t)
Definition: SSEVec.h:18
double pt() const
track transverse momentum
Definition: TrackBase.h:660
double chi2() const
chi-squares
Definition: Vertex.h:98
double z() const
z coordinate
Definition: Vertex.h:115
math::XYZPoint Point
point in the space
Definition: Vertex.h:39
float trackWeight(const TREF &r) const
returns the weight with which a Track has contributed to the vertex-fit.
Definition: Vertex.h:81
math::Error< dimension4D >::type CovarianceMatrix4D
covariance error matrix (4x4)
Definition: Vertex.h:49
double ndof() const
Definition: Vertex.h:105
CovarianceMatrix covariance() const
return SMatrix
Definition: Vertex.h:134
void add(const TrackBaseRef &r, float w=1.0)
add a reference to a Track
double x() const
x coordinate
Definition: Vertex.h:111
double xError() const
error on x
Definition: Vertex.h:119
bool isFake() const
Definition: Vertex.h:72
unsigned int nTracks(float minWeight=0.5) const
Returns the number of tracks in the vertex with weight above minWeight.
Definition: Vertex.cc:169
Point position_
position
Definition: Vertex.h:186
Error4D error4D() const
return SMatrix
Definition: Vertex.h:141
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
math::XYZTLorentzVectorD p4(float mass=0.13957018, float minWeight=0.5) const
Returns the four momentum of the sum of the tracks, assuming the given mass for the decay products...
Definition: Vertex.cc:135
double b
Definition: hdecay.h:120
Error error() const
return SMatrix
Definition: Vertex.h:139
bool operator()(const Track &t) const
Definition: Vertex.h:175
float covariance_[size4D]
covariance matrix (4x4) as vector
Definition: Vertex.h:188
fixed size matrix
double a
Definition: hdecay.h:121
double time_
Definition: Vertex.h:196
math::Error< dimension4D >::type Error4D
covariance error matrix (4x4)
Definition: Vertex.h:47
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector<TrackRef>
Definition: Vertex.h:37
double normalizedChi2() const
chi-squared divided by n.d.o.f.
Definition: Vertex.h:107
void removeTracks()
Definition: Vertex.cc:101
trackRef_iterator tracks_begin() const
first iterator over tracks
Definition: Vertex.cc:76
index idx(index i, index j) const
position index
Definition: Vertex.h:199
double yError() const
error on y
Definition: Vertex.h:121
double tError() const
error on t
Definition: Vertex.h:125
size_t tracksSize() const
number of tracks
Definition: Vertex.cc:71
double t() const
t coordinate
Definition: Vertex.h:117
float ndof_
number of degrees of freedom
Definition: Vertex.h:184
bool validity_
tells wether the vertex is really valid.
Definition: Vertex.h:195
std::vector< uint8_t > weights_
Definition: Vertex.h:193