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>
27 #include <Math/GenVector/PxPyPzE4D.h>
28 #include <Math/GenVector/PxPyPzM4D.h>
30 
31 namespace reco {
32 
33  class Track;
34 
35  class Vertex {
36  public:
38  typedef std::vector<TrackBaseRef>::const_iterator trackRef_iterator;
42  enum { dimension = 3, dimension4D = 4 };
52  enum { size = dimension * (dimension + 1) / 2, size4D = (dimension4D) * (dimension4D + 1) / 2 };
54  typedef unsigned int index;
58  Vertex() : chi2_(0.0), ndof_(0), position_(0., 0., 0.), time_(0.) {
59  validity_ = false;
60  for (int i = 0; i < size4D; ++i)
61  covariance_[i] = 0.;
62  }
64  Vertex(const Point &, const Error &);
66  Vertex(const Point &, const Error4D &, double);
68  Vertex(const Point &, const Error &, double chi2, double ndof, size_t size);
70  Vertex(const Point &, const Error4D &, double time, double chi2, double ndof, size_t size);
72  bool isValid() const { return validity_; }
76  bool isFake() const { return (chi2_ == 0 && ndof_ == 0 && tracks_.empty()); }
78  void add(const TrackBaseRef &r, float w = 1.0);
80  void add(const TrackBaseRef &r, const Track &refTrack, float w = 1.0);
81  void removeTracks();
82 
84  template <typename TREF>
85  float trackWeight(const TREF &r) const {
86  int i = 0;
87  for (auto const &t : tracks_) {
88  if ((r.id() == t.id()) & (t.key() == r.key()))
89  return weights_[i] / 255.f;
90  ++i;
91  }
92  return 0;
93  }
99  size_t tracksSize() const;
101  const TrackBaseRef &trackRefAt(size_t idx) const { return tracks_[idx]; }
103  double chi2() const { return chi2_; }
110  double ndof() const { return ndof_; }
112  double normalizedChi2() const { return ndof_ != 0 ? chi2_ / ndof_ : chi2_ * 1e6; }
114  const Point &position() const { return position_; }
116  double x() const { return position_.X(); }
118  double y() const { return position_.Y(); }
120  double z() const { return position_.Z(); }
122  double t() const { return time_; }
124  double xError() const { return sqrt(covariance(0, 0)); }
126  double yError() const { return sqrt(covariance(1, 1)); }
128  double zError() const { return sqrt(covariance(2, 2)); }
130  double tError() const { return sqrt(covariance(3, 3)); }
132  // Note that:
133  // double error( int i, int j ) const
134  // is OBSOLETE, use covariance(i, j)
135  double covariance(int i, int j) const { return covariance_[idx(i, j)]; }
138  Error m;
139  fill(m);
140  return m;
141  }
144  Error4D m;
145  fill(m);
146  return m;
147  }
148 
150  Error error() const {
151  Error m;
152  fill(m);
153  return m;
154  }
156  Error4D error4D() const {
157  Error4D m;
158  fill(m);
159  return m;
160  }
161 
163  void fill(CovarianceMatrix &v) const;
165  void fill(CovarianceMatrix4D &v) const;
166 
168  bool hasRefittedTracks() const { return !refittedTracks_.empty(); }
169 
172  TrackBaseRef originalTrack(const Track &refTrack) const;
173 
176  Track refittedTrack(const TrackBaseRef &track) const;
177 
180  Track refittedTrack(const TrackRef &track) const;
181 
183  const std::vector<Track> &refittedTracks() const { return refittedTracks_; }
184 
186  math::XYZTLorentzVectorD p4(float mass = 0.13957018, float minWeight = 0.5) const;
187 
189  unsigned int nTracks(float minWeight = 0.5) const;
190 
191  class TrackEqual {
192  public:
193  TrackEqual(const Track &t) : track_(t) {}
194  bool operator()(const Track &t) const { return t.pt() == track_.pt(); }
195 
196  private:
197  const Track &track_;
198  };
199 
200  private:
202  float chi2_;
204  float ndof_;
210  std::vector<TrackBaseRef> tracks_;
212  std::vector<Track> refittedTracks_;
213  std::vector<uint8_t> weights_;
215  bool validity_;
216  double time_;
217 
219  index idx(index i, index j) const {
220  int a = (i <= j ? i : j), b = (i <= j ? j : i);
221  return b * (b + 1) / 2 + a;
222  }
223  };
224 
225 } // namespace reco
226 
227 #endif
reco::Vertex::trackRef_iterator
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector<TrackRef>
Definition: Vertex.h:38
reco::Vertex::isValid
bool isValid() const
Tells whether the vertex is valid.
Definition: Vertex.h:72
reco::Vertex::time_
double time_
Definition: Vertex.h:216
reco::Vertex::originalTrack
TrackBaseRef originalTrack(const Track &refTrack) const
Definition: Vertex.cc:96
reco::Vertex::TrackEqual::track_
const Track & track_
Definition: Vertex.h:197
mps_fire.i
i
Definition: mps_fire.py:355
reco::Vertex::chi2_
float chi2_
chi-sqared
Definition: Vertex.h:202
reco::Vertex::z
double z() const
z coordinate
Definition: Vertex.h:120
reco::Vertex::Error
math::Error< dimension >::type Error
covariance error matrix (3x3)
Definition: Vertex.h:44
math::XYZTLorentzVectorD
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< double > > XYZTLorentzVectorD
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:14
reco::Vertex::refittedTrack
Track refittedTrack(const TrackBaseRef &track) const
Definition: Vertex.cc:106
reco::Vertex::position
const Point & position() const
position
Definition: Vertex.h:114
CovarianceMatrix
Definition: CovarianceMatrix.h:27
reco::Vertex::idx
index idx(index i, index j) const
position index
Definition: Vertex.h:219
reco::Vertex::CovarianceMatrix
math::Error< dimension >::type CovarianceMatrix
covariance error matrix (3x3)
Definition: Vertex.h:46
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:45
reco::Vertex::validity_
bool validity_
tells wether the vertex is really valid.
Definition: Vertex.h:215
findQualityFiles.v
v
Definition: findQualityFiles.py:179
fftjetcommon_cfi.minWeight
minWeight
Definition: fftjetcommon_cfi.py:213
reco::Vertex::TrackEqual
Definition: Vertex.h:191
reco::Vertex::tError
double tError() const
error on t
Definition: Vertex.h:130
edm::Ref< TrackCollection >
reco::TrackBase::pt
double pt() const
track transverse momentum
Definition: TrackBase.h:608
reco::Vertex::dimension4D
Definition: Vertex.h:42
reco::Vertex::tracks_
std::vector< TrackBaseRef > tracks_
reference to tracks
Definition: Vertex.h:210
Track.h
reco::Vertex::isFake
bool isFake() const
Definition: Vertex.h:76
TrackFwd.h
reco::Vertex::xError
double xError() const
error on x
Definition: Vertex.h:124
w
const double w
Definition: UKUtility.cc:23
visualization-live-secondInstance_cfg.m
m
Definition: visualization-live-secondInstance_cfg.py:72
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
reco::Vertex::tracks_end
trackRef_iterator tracks_end() const
last iterator over tracks
Definition: Vertex.cc:73
reco::Vertex::add
void add(const TrackBaseRef &r, float w=1.0)
add a reference to a Track
Definition: Vertex.cc:79
reco::Track
Definition: Track.h:27
b
double b
Definition: hdecay.h:118
RefToBase.h
reco::Vertex::refittedTracks
const std::vector< Track > & refittedTracks() const
Returns the container of refitted tracks.
Definition: Vertex.h:183
reco::Vertex::hasRefittedTracks
bool hasRefittedTracks() const
Checks whether refitted tracks are stored.
Definition: Vertex.h:168
reco::Vertex::zError
double zError() const
error on z
Definition: Vertex.h:128
reco::Vertex::covariance
CovarianceMatrix covariance() const
return SMatrix
Definition: Vertex.h:137
Error.h
Point
Structure Point Contains parameters of Gaussian fits to DMRs.
Definition: DMRtrends.cc:57
reco::Vertex::p4
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
reco::Vertex::refittedTracks_
std::vector< Track > refittedTracks_
The vector of refitted tracks.
Definition: Vertex.h:212
math::XYZPoint
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
a
double a
Definition: hdecay.h:119
reco::Vertex::dimension
Definition: Vertex.h:42
reco::Vertex::tracks_begin
trackRef_iterator tracks_begin() const
first iterator over tracks
Definition: Vertex.cc:71
reco::Vertex::position_
Point position_
position
Definition: Vertex.h:206
reco::Vertex::x
double x() const
x coordinate
Definition: Vertex.h:116
reco::Vertex::CovarianceMatrix4D
math::Error< dimension4D >::type CovarianceMatrix4D
covariance error matrix (4x4)
Definition: Vertex.h:50
LorentzVector.h
reco::Vertex::covariance4D
CovarianceMatrix4D covariance4D() const
return SMatrix 4D
Definition: Vertex.h:143
reco::Vertex::tracksSize
size_t tracksSize() const
number of tracks
Definition: Vertex.cc:69
reco::Vertex::error
Error error() const
return SMatrix
Definition: Vertex.h:150
reco::Vertex::chi2
double chi2() const
chi-squares
Definition: Vertex.h:103
alignCSCRings.r
r
Definition: alignCSCRings.py:93
reco::Vertex::TrackEqual::operator()
bool operator()(const Track &t) const
Definition: Vertex.h:194
VertexFwd.h
reco::Vertex::trackWeight
float trackWeight(const TREF &r) const
returns the weight with which a Track has contributed to the vertex-fit.
Definition: Vertex.h:85
reco::Vertex::size
Definition: Vertex.h:52
reco::Vertex::TrackEqual::TrackEqual
TrackEqual(const Track &t)
Definition: Vertex.h:193
reco::Vertex::covariance
double covariance(int i, int j) const
(i, j)-th element of error matrix, i, j = 0, ... 2
Definition: Vertex.h:135
reco::Vertex::Point
math::XYZPoint Point
point in the space
Definition: Vertex.h:40
reco::Vertex::y
double y() const
y coordinate
Definition: Vertex.h:118
reco::Vertex::error4D
Error4D error4D() const
return SMatrix
Definition: Vertex.h:156
reco::Vertex::fill
void fill(CovarianceMatrix &v) const
fill SMatrix
reco::Vertex::t
double t() const
t coordinate
Definition: Vertex.h:122
math::Error::type
ErrorD< N >::type type
Definition: Error.h:32
reco::Vertex::normalizedChi2
double normalizedChi2() const
chi-squared divided by n.d.o.f.
Definition: Vertex.h:112
reco::Vertex::Vertex
Vertex()
Definition: Vertex.h:58
reco::Vertex::covariance_
float covariance_[size4D]
covariance matrix (4x4) as vector
Definition: Vertex.h:208
EgHLTOffHistBins_cfi.mass
mass
Definition: EgHLTOffHistBins_cfi.py:34
Point3D.h
edm::RefToBase< reco::Track >
reco::Vertex::ndof_
float ndof_
number of degrees of freedom
Definition: Vertex.h:204
HLT_2018_cff.track
track
Definition: HLT_2018_cff.py:10352
reco::Vertex::nTracks
unsigned int nTracks(float minWeight=0.5) const
Returns the number of tracks in the vertex with weight above minWeight.
Definition: Vertex.cc:146
reco::Vertex::yError
double yError() const
error on y
Definition: Vertex.h:126
reco::Vertex::removeTracks
void removeTracks()
Definition: Vertex.cc:90
reco::Vertex::index
unsigned int index
index type
Definition: Vertex.h:54
dqmiolumiharvest.j
j
Definition: dqmiolumiharvest.py:66
ntuplemaker.time
time
Definition: ntuplemaker.py:310
reco::Vertex::ndof
double ndof() const
Definition: Vertex.h:110
reco::Vertex
Definition: Vertex.h:35
reco::Vertex::trackRefAt
const TrackBaseRef & trackRefAt(size_t idx) const
python friendly track getting
Definition: Vertex.h:101
reco::Vertex::Error4D
math::Error< dimension4D >::type Error4D
covariance error matrix (4x4)
Definition: Vertex.h:48
reco::Vertex::weights_
std::vector< uint8_t > weights_
Definition: Vertex.h:213
math::Error
fixed size error matrix
Definition: Error.h:31
reco::Vertex::size4D
Definition: Vertex.h:52