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, 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.) {
58  validity_ = false;
59  for (int i = 0; i < size4D; ++i)
60  covariance_[i] = 0.;
61  }
63  Vertex(const Point &, const Error &);
65  Vertex(const Point &, const Error4D &, double);
67  Vertex(const Point &, const Error &, double chi2, double ndof, size_t size);
69  Vertex(const Point &, const Error4D &, double time, double chi2, double ndof, size_t size);
71  bool isValid() const { return validity_; }
75  bool isFake() const { return (chi2_ == 0 && ndof_ == 0 && tracks_.empty()); }
77  void add(const TrackBaseRef &r, float w = 1.0);
79  void add(const TrackBaseRef &r, const Track &refTrack, float w = 1.0);
80  void removeTracks();
81 
83  template <typename TREF>
84  float trackWeight(const TREF &r) const {
85  int i = 0;
86  for (auto const &t : tracks_) {
87  if ((r.id() == t.id()) & (t.key() == r.key()))
88  return weights_[i] / 255.f;
89  ++i;
90  }
91  return 0;
92  }
94  trackRef_iterator tracks_begin() const;
96  trackRef_iterator tracks_end() const;
98  size_t tracksSize() const;
100  const TrackBaseRef &trackRefAt(size_t idx) const { return tracks_[idx]; }
102  double chi2() const { return chi2_; }
109  double ndof() const { return ndof_; }
111  double normalizedChi2() const { return ndof_ != 0 ? chi2_ / ndof_ : chi2_ * 1e6; }
113  const Point &position() const { return position_; }
115  double x() const { return position_.X(); }
117  double y() const { return position_.Y(); }
119  double z() const { return position_.Z(); }
121  double t() const { return time_; }
123  double xError() const { return sqrt(covariance(0, 0)); }
125  double yError() const { return sqrt(covariance(1, 1)); }
127  double zError() const { return sqrt(covariance(2, 2)); }
129  double tError() const { return sqrt(covariance(3, 3)); }
131  // Note that:
132  // double error( int i, int j ) const
133  // is OBSOLETE, use covariance(i, j)
134  double covariance(int i, int j) const { return covariance_[idx(i, j)]; }
136  CovarianceMatrix covariance() const {
137  Error m;
138  fill(m);
139  return m;
140  }
142  CovarianceMatrix4D covariance4D() const {
143  Error4D m;
144  fill(m);
145  return m;
146  }
147 
149  Error error() const {
150  Error m;
151  fill(m);
152  return m;
153  }
155  Error4D error4D() const {
156  Error4D m;
157  fill(m);
158  return m;
159  }
160 
162  void fill(CovarianceMatrix &v) const;
164  void fill(CovarianceMatrix4D &v) const;
165 
167  bool hasRefittedTracks() const { return !refittedTracks_.empty(); }
168 
171  TrackBaseRef originalTrack(const Track &refTrack) const;
172 
175  Track refittedTrack(const TrackBaseRef &track) const;
176 
179  Track refittedTrack(const TrackRef &track) const;
180 
182  const std::vector<Track> &refittedTracks() const { return refittedTracks_; }
183 
185  math::XYZTLorentzVectorD p4(float mass = 0.13957018, float minWeight = 0.5) const;
186 
188  unsigned int nTracks(float minWeight = 0.5) const;
189 
190  class TrackEqual {
191  public:
192  TrackEqual(const Track &t) : track_(t) {}
193  bool operator()(const Track &t) const { return t.pt() == track_.pt(); }
194 
195  private:
196  const Track &track_;
197  };
198 
199  private:
201  float chi2_;
203  float ndof_;
205  Point position_;
209  std::vector<TrackBaseRef> tracks_;
211  std::vector<Track> refittedTracks_;
212  std::vector<uint8_t> weights_;
214  bool validity_;
215  double time_;
216 
218  index idx(index i, index j) const {
219  int a = (i <= j ? i : j), b = (i <= j ? j : i);
220  return b * (b + 1) / 2 + a;
221  }
222  };
223 
224 } // namespace reco
225 
226 #endif
void fill(CovarianceMatrix &v) const
fill SMatrix
float chi2_
chi-sqared
Definition: Vertex.h:201
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:73
double zError() const
error on z
Definition: Vertex.h:127
TrackBaseRef originalTrack(const Track &refTrack) const
Definition: Vertex.cc:96
double y() const
y coordinate
Definition: Vertex.h:117
bool isValid() const
Tells whether the vertex is valid.
Definition: Vertex.h:71
std::vector< Track > refittedTracks_
The vector of refitted tracks.
Definition: Vertex.h:211
bool hasRefittedTracks() const
Checks whether refitted tracks are stored.
Definition: Vertex.h:167
const TrackBaseRef & trackRefAt(size_t idx) const
python friendly track getting
Definition: Vertex.h:100
math::Error< dimension >::type Error
covariance error matrix (3x3)
Definition: Vertex.h:43
std::vector< TrackBaseRef > tracks_
reference to tracks
Definition: Vertex.h:209
ErrorD< N >::type type
Definition: Error.h:32
CovarianceMatrix4D covariance4D() const
return SMatrix 4D
Definition: Vertex.h:142
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:134
const std::vector< Track > & refittedTracks() const
Returns the container of refitted tracks.
Definition: Vertex.h:182
const Point & position() const
position
Definition: Vertex.h:113
math::Error< dimension >::type CovarianceMatrix
covariance error matrix (3x3)
Definition: Vertex.h:45
TrackEqual(const Track &t)
Definition: Vertex.h:192
const Track & track_
Definition: Vertex.h:196
T sqrt(T t)
Definition: SSEVec.h:19
double pt() const
track transverse momentum
Definition: TrackBase.h:602
double chi2() const
chi-squares
Definition: Vertex.h:102
double z() const
z coordinate
Definition: Vertex.h:119
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:84
math::Error< dimension4D >::type CovarianceMatrix4D
covariance error matrix (4x4)
Definition: Vertex.h:49
void add(const TrackBaseRef &r, float w=1.0)
add a reference to a Track
Definition: Vertex.cc:79
double ndof() const
Definition: Vertex.h:109
CovarianceMatrix covariance() const
return SMatrix
Definition: Vertex.h:136
double x() const
x coordinate
Definition: Vertex.h:115
double xError() const
error on x
Definition: Vertex.h:123
bool isFake() const
Definition: Vertex.h:75
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:205
Error4D error4D() const
return SMatrix
Definition: Vertex.h:155
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:118
double b
Definition: hdecay.h:118
Error error() const
return SMatrix
Definition: Vertex.h:149
bool operator()(const Track &t) const
Definition: Vertex.h:193
float covariance_[size4D]
covariance matrix (4x4) as vector
Definition: Vertex.h:207
fixed size matrix
Track refittedTrack(const TrackBaseRef &track) const
Definition: Vertex.cc:106
double a
Definition: hdecay.h:119
double time_
Definition: Vertex.h:215
math::Error< dimension4D >::type Error4D
covariance error matrix (4x4)
Definition: Vertex.h:47
double normalizedChi2() const
chi-squared divided by n.d.o.f.
Definition: Vertex.h:111
void removeTracks()
Definition: Vertex.cc:90
trackRef_iterator tracks_begin() const
first iterator over tracks
Definition: Vertex.cc:71
index idx(index i, index j) const
position index
Definition: Vertex.h:218
double yError() const
error on y
Definition: Vertex.h:125
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector<TrackRef>
Definition: Vertex.h:37
double tError() const
error on t
Definition: Vertex.h:129
size_t tracksSize() const
number of tracks
Definition: Vertex.cc:69
double t() const
t coordinate
Definition: Vertex.h:121
float ndof_
number of degrees of freedom
Definition: Vertex.h:203
bool validity_
tells wether the vertex is really valid.
Definition: Vertex.h:214
std::vector< uint8_t > weights_
Definition: Vertex.h:212