CMS 3D CMS Logo

All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 };
47  enum { size = dimension * ( dimension + 1 ) / 2 };
49  typedef unsigned int index;
53  Vertex(): chi2_( 0.0 ), ndof_( 0 ), position_(0.,0.,0. ) { validity_ = false; for(int i = 0; i < size; ++ i ) covariance_[i]=0.;
54 }
56  Vertex( const Point &, const Error &);
58  Vertex( const Point &, const Error &, double chi2, double ndof, size_t size );
60  bool isValid() const {return validity_;}
64  bool isFake() const {return (chi2_==0 && ndof_==0 && tracks_.empty());}
66  void add( const TrackBaseRef & r, float w=1.0 );
68  void add( const TrackBaseRef & r, const Track & refTrack, float w=1.0 );
69  void removeTracks();
71  float trackWeight ( const TrackBaseRef & r ) const;
73  float trackWeight ( const TrackRef & r ) const;
79  size_t tracksSize() const;
81  double chi2() const { return chi2_; }
88  double ndof() const { return ndof_; }
90  double normalizedChi2() const { return ndof_ != 0 ? chi2_ / ndof_ : chi2_ * 1e6; }
92  const Point & position() const { return position_; }
94  double x() const { return position_.X(); }
96  double y() const { return position_.Y(); }
98  double z() const { return position_.Z(); }
100  double xError() const { return sqrt( covariance(0, 0) ); }
102  double yError() const { return sqrt( covariance(1, 1) ); }
104  double zError() const { return sqrt( covariance(2, 2) ); }
106  // Note that:
107  // double error( int i, int j ) const
108  // is OBSOLETE, use covariance(i, j)
109  double covariance( int i, int j ) const {
110  return covariance_[ idx( i, j ) ];
111  }
113  CovarianceMatrix covariance() const { Error m; fill( m ); return m; }
115  Error error() const { Error m; fill( m ); return m; }
117  void fill( CovarianceMatrix & v ) const;
118 
120  bool hasRefittedTracks() const {return !refittedTracks_.empty();}
121 
124  TrackBaseRef originalTrack(const Track & refTrack) const;
125 
128  Track refittedTrack(const TrackBaseRef & track) const;
129 
132  Track refittedTrack(const TrackRef & track) const;
133 
135  const std::vector<Track> & refittedTracks() const { return refittedTracks_;}
136 
138  math::XYZTLorentzVectorD p4(float mass=0.13957018,float minWeight=0.5) const;
139 
141  unsigned int nTracks(float minWeight=0.5) const;
142 
143  class TrackEqual {
144  public:
145  TrackEqual( const Track & t) : track_( t ) { }
146  bool operator()( const Track & t ) const { return t.pt()==track_.pt();}
147  private:
148  const Track & track_;
149  };
150 
151  private:
153  float chi2_;
155  float ndof_;
159  float covariance_[ size ];
161  std::vector<TrackBaseRef > tracks_;
163  std::vector<Track> refittedTracks_;
164  std::vector<uint8_t> weights_;
166  bool validity_;
167 
168 
170  index idx( index i, index j ) const {
171  int a = ( i <= j ? i : j ), b = ( i <= j ? j : i );
172  return b * ( b + 1 ) / 2 + a;
173  }
174  };
175 
176 }
177 
178 #endif
int i
Definition: DBlmapReader.cc:9
std::vector< TrackBaseRef > tracks_
reference to tracks
Definition: Vertex.h:161
float chi2_
chi-sqared
Definition: Vertex.h:153
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< double > > XYZTLorentzVectorD
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:14
trackRef_iterator tracks_end() const
last iterator over tracks
Definition: Vertex.cc:44
double zError() const
error on z
Definition: Vertex.h:104
TrackBaseRef originalTrack(const Track &refTrack) const
Definition: Vertex.cc:85
Track refittedTrack(const TrackBaseRef &track) const
double y() const
y coordinate
Definition: Vertex.h:96
bool isValid() const
Tells whether the vertex is valid.
Definition: Vertex.h:60
std::vector< Track > refittedTracks_
The vector of refitted tracks.
Definition: Vertex.h:163
bool hasRefittedTracks() const
Checks whether refitted tracks are stored.
Definition: Vertex.h:120
fixed size error matrix
Definition: Error.h:28
math::Error< dimension >::type Error
covariance error matrix (3x3)
Definition: Vertex.h:43
void fill(CovarianceMatrix &v) const
fill SMatrix
Definition: Vertex.cc:27
ErrorD< N >::type type
Definition: Error.h:29
unsigned int index
index type
Definition: Vertex.h:49
double covariance(int i, int j) const
(i, j)-th element of error matrix, i, j = 0, ... 2
Definition: Vertex.h:109
const std::vector< Track > & refittedTracks() const
Returns the container of refitted tracks.
Definition: Vertex.h:135
const Point & position() const
position
Definition: Vertex.h:92
math::Error< dimension >::type CovarianceMatrix
covariance error matrix (3x3)
Definition: Vertex.h:45
TrackEqual(const Track &t)
Definition: Vertex.h:145
const Track & track_
Definition: Vertex.h:148
float trackWeight(const TrackBaseRef &r) const
returns the weight with which a Track has contributed to the vertex-fit.
T sqrt(T t)
Definition: SSEVec.h:48
double pt() const
track transverse momentum
Definition: TrackBase.h:129
math::XYZPoint Point
double chi2() const
chi-squares
Definition: Vertex.h:81
int j
Definition: DBlmapReader.cc:9
double z() const
y coordinate
Definition: Vertex.h:98
math::XYZPoint Point
point in the space
Definition: Vertex.h:39
double ndof() const
Definition: Vertex.h:88
CovarianceMatrix covariance() const
return SMatrix
Definition: Vertex.h:113
void add(const TrackBaseRef &r, float w=1.0)
add a reference to a Track
double x() const
x coordinate
Definition: Vertex.h:94
double xError() const
error on x
Definition: Vertex.h:100
bool isFake() const
Definition: Vertex.h:64
unsigned int nTracks(float minWeight=0.5) const
Returns the number of tracks in the vertex with weight above minWeight.
Definition: Vertex.cc:146
Point position_
position
Definition: Vertex.h:157
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:112
double b
Definition: hdecay.h:120
Error error() const
return SMatrix
Definition: Vertex.h:115
bool operator()(const Track &t) const
Definition: Vertex.h:146
double a
Definition: hdecay.h:121
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector&lt;TrackRef&gt;
Definition: Vertex.h:37
double normalizedChi2() const
chi-squared divided by n.d.o.f.
Definition: Vertex.h:90
T w() const
void removeTracks()
Definition: Vertex.cc:64
trackRef_iterator tracks_begin() const
first iterator over tracks
Definition: Vertex.cc:39
float covariance_[size]
covariance matrix (3x3) as vector
Definition: Vertex.h:159
index idx(index i, index j) const
position index
Definition: Vertex.h:170
double yError() const
error on y
Definition: Vertex.h:102
size_t tracksSize() const
number of tracks
Definition: Vertex.cc:34
float ndof_
number of degrees of freedom
Definition: Vertex.h:155
bool validity_
tells wether the vertex is really valid.
Definition: Vertex.h:166
std::vector< uint8_t > weights_
Definition: Vertex.h:164