CMS 3D CMS Logo

Vertex.cc
Go to the documentation of this file.
2 #include <Math/GenVector/PxPyPzE4D.h>
3 #include <Math/GenVector/PxPyPzM4D.h>
4 
5 using namespace reco;
6 using namespace std;
7 
8 Vertex::Vertex( const Point & p , const Error & err, double chi2, double ndof, size_t size ) :
9  chi2_( chi2 ), ndof_( ndof ), position_( p ), time_(0.) {
10  tracks_.reserve( size );
11  index idx = 0;
12  for( index i = 0; i < dimension4D; ++ i ) {
13  for( index j = 0; j <= i; ++ j ) {
14  if( i == dimension || j == dimension ) {
15  covariance_[ idx ++ ] = 0.0;
16  } else {
17  covariance_[ idx ++ ] = err( i, j );
18  }
19  }
20  }
21  validity_ = true;
22 }
23 
24 Vertex::Vertex( const Point & p , const Error4D & err, double time, double chi2, double ndof, size_t size ) :
25  chi2_( chi2 ), ndof_( ndof ), position_( p ), time_(time) {
26  tracks_.reserve( size4D );
27  index idx = 0;
28  for( index i = 0; i < dimension4D; ++ i )
29  for( index j = 0; j <= i; ++ j )
30  covariance_[ idx ++ ] = err( i, j );
31  validity_ = true;
32 }
33 
34 Vertex::Vertex( const Point & p , const Error & err) :
35  chi2_( 0.0 ), ndof_( 0 ), position_( p ), time_(0.) {
36  index idx = 0;
37  for( index i = 0; i < dimension4D; ++ i ) {
38  for( index j = 0; j <= i; ++ j ) {
39  if( i == dimension || j == dimension ) {
40  covariance_[ idx ++ ] = 0.0;
41  } else {
42  covariance_[ idx ++ ] = err( i, j );
43  }
44  }
45  }
46  validity_ = true;
47 }
48 
49 Vertex::Vertex( const Point & p , const Error4D & err, double time) :
50  chi2_( 0.0 ), ndof_( 0 ), position_( p ), time_(time) {
51  index idx = 0;
52  for( index i = 0; i < dimension + 1; ++ i )
53  for( index j = 0; j <= i; ++ j )
54  covariance_[ idx ++ ] = err( i, j );
55  validity_ = true;
56 }
57 
58 void Vertex::fill( Error & err ) const {
59  Error4D temp;
60  fill(temp);
61  err = temp.Sub<Error>(0,0);
62 }
63 
64 void Vertex::fill( Error4D & err ) const {
65  index idx = 0;
66  for( index i = 0; i < dimension4D; ++ i )
67  for( index j = 0; j <= i; ++ j )
68  err( i, j ) = covariance_[ idx ++ ];
69 }
70 
71 size_t Vertex::tracksSize() const
72 {
73  return weights_.size();
74 }
75 
77 {
78  return tracks_.begin();
79 }
80 
82 {
83 // if ( !(tracks_.size() ) ) createTracks();
84  return tracks_.end();
85  // return weights_.keys().end();
86 }
87 
88 void Vertex::add ( const TrackBaseRef & r, float w )
89 {
90  tracks_.push_back ( r );
91  weights_.push_back(w*255);
92 }
93 
94 void Vertex::add ( const TrackBaseRef & r, const Track & refTrack, float w )
95 {
96  tracks_.push_back ( r );
97  refittedTracks_.push_back ( refTrack );
98  weights_.push_back(w*255);
99 }
100 
102 {
103  weights_.clear();
104  tracks_.clear();
105  refittedTracks_.clear();
106 }
107 
108 TrackBaseRef Vertex::originalTrack(const Track & refTrack) const
109 {
110  if (refittedTracks_.empty())
111  throw cms::Exception("Vertex") << "No refitted tracks stored in vertex\n";
112  std::vector<Track>::const_iterator it =
113  find_if(refittedTracks_.begin(), refittedTracks_.end(), TrackEqual(refTrack));
114  if (it==refittedTracks_.end())
115  throw cms::Exception("Vertex") << "Refitted track not found in list\n";
116  size_t pos = it - refittedTracks_.begin();
117  return tracks_[pos];
118 }
119 
121 {
122  if (refittedTracks_.empty())
123  throw cms::Exception("Vertex") << "No refitted tracks stored in vertex\n";
124  trackRef_iterator it = find(tracks_begin(), tracks_end(), track);
125  if (it==tracks_end()) throw cms::Exception("Vertex") << "Track not found in list\n";
126  size_t pos = it - tracks_begin();
127  return refittedTracks_[pos];
128 }
129 
131 {
132  return refittedTrack(TrackBaseRef(track));
133 }
134 
136 {
137 
139  ROOT::Math::LorentzVector<ROOT::Math::PxPyPzM4D<double> > vec;
140 
141  if(hasRefittedTracks()) {
142  for(std::vector<Track>::const_iterator iter = refittedTracks_.begin();
143  iter != refittedTracks_.end(); ++iter) {
144  if (trackWeight(originalTrack(*iter)) >=minWeight) {
145  vec.SetPx(iter->px());
146  vec.SetPy(iter->py());
147  vec.SetPz(iter->pz());
148  vec.SetM(mass);
149  sum += vec;
150  }
151  }
152  }
153  else
154  {
155  for(std::vector<reco::TrackBaseRef>::const_iterator iter = tracks_begin();
156  iter != tracks_end(); iter++) {
157  if (trackWeight(*iter) >=minWeight) {
158  vec.SetPx((*iter)->px());
159  vec.SetPy((*iter)->py());
160  vec.SetPz((*iter)->pz());
161  vec.SetM(mass);
162  sum += vec;
163  }
164  }
165  }
166  return sum;
167 }
168 
169 unsigned int Vertex::nTracks(float minWeight) const
170 {
171  int n=0;
172  if(hasRefittedTracks()) {
173  for(std::vector<Track>::const_iterator iter = refittedTracks_.begin(); iter != refittedTracks_.end(); ++iter)
174  if (trackWeight(originalTrack(*iter)) >=minWeight)
175  n++;
176  }
177  else
178  {
179  for(std::vector<reco::TrackBaseRef>::const_iterator iter = tracks_begin(); iter != tracks_end(); iter++)
180  if (trackWeight(*iter) >=minWeight)
181  n++;
182  }
183  return n;
184 }
185 
size
Write out results.
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
TrackBaseRef originalTrack(const Track &refTrack) const
Definition: Vertex.cc:108
Track refittedTrack(const TrackBaseRef &track) const
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
math::Error< dimension >::type Error
covariance error matrix (3x3)
Definition: Vertex.h:43
std::vector< TrackBaseRef > tracks_
reference to tracks
Definition: Vertex.h:190
unsigned int index
index type
Definition: Vertex.h:53
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
edm::RefToBase< reco::Track > TrackBaseRef
persistent reference to a Track, using views
Definition: TrackFwd.h:36
double chi2() const
chi-squares
Definition: Vertex.h:98
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
double ndof() const
Definition: Vertex.h:105
void add(const TrackBaseRef &r, float w=1.0)
add a reference to a Track
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
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
float covariance_[size4D]
covariance matrix (4x4) as vector
Definition: Vertex.h:188
fixed size matrix
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
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
size_t tracksSize() const
number of tracks
Definition: Vertex.cc:71
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