CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
Vertex.cc
Go to the documentation of this file.
2 #include <Math/GenVector/PxPyPzE4D.h>
3 #include <Math/GenVector/PxPyPzM4D.h>
4 
5 // $Id: Vertex.cc,v 1.17 2010/04/16 08:08:28 arizzi Exp $
6 using namespace reco;
7 using namespace std;
8 
9 Vertex::Vertex( const Point & p , const Error & err, double chi2, double ndof, size_t size ) :
10  chi2_( chi2 ), ndof_( ndof ), position_( p ) {
11  tracks_.reserve( size );
12  index idx = 0;
13  for( index i = 0; i < dimension; ++ i )
14  for( index j = 0; j <= i; ++ j )
15  covariance_[ idx ++ ] = err( i, j );
16  validity_ = true;
17 }
18 
19 Vertex::Vertex( const Point & p , const Error & err) :
20  chi2_( 0.0 ), ndof_( 0 ), position_( p ) {
21  index idx = 0;
22  for( index i = 0; i < dimension; ++ i )
23  for( index j = 0; j <= i; ++ j )
24  covariance_[ idx ++ ] = err( i, j );
25  validity_ = true;
26 }
27 
28 void Vertex::fill( Error & err ) const {
29  index idx = 0;
30  for( index i = 0; i < dimension; ++ i )
31  for( index j = 0; j <= i; ++ j )
32  err( i, j ) = covariance_[ idx ++ ];
33 }
34 
35 size_t Vertex::tracksSize() const
36 {
37  return weights_.size();
38 }
39 
41 {
42  return tracks_.begin();
43 }
44 
46 {
47 // if ( !(tracks_.size() ) ) createTracks();
48  return tracks_.end();
49  // return weights_.keys().end();
50 }
51 
52 void Vertex::add ( const TrackBaseRef & r, float w )
53 {
54  tracks_.push_back ( r );
55  weights_.push_back(w);
56 
57 }
58 
59 void Vertex::add ( const TrackBaseRef & r, const Track & refTrack, float w )
60 {
61  tracks_.push_back ( r );
62  refittedTracks_.push_back ( refTrack );
63  weights_.push_back(w);
64 }
65 
67 {
68  weights_.clear();
69  tracks_.clear();
70  refittedTracks_.clear();
71 }
72 
73 float Vertex::trackWeight ( const TrackBaseRef & track ) const
74 {
76  if (it==tracks_end()) return 0.0;
77  size_t pos = it - tracks_begin();
78  return weights_[pos];
79 }
80 
81 float Vertex::trackWeight ( const TrackRef & track ) const
82 {
83  return trackWeight(TrackBaseRef(track));
84 }
85 
86 
87 TrackBaseRef Vertex::originalTrack(const Track & refTrack) const
88 {
89  if (refittedTracks_.empty())
90  throw cms::Exception("Vertex") << "No refitted tracks stored in vertex\n";
91  std::vector<Track>::const_iterator it =
92  find_if(refittedTracks_.begin(), refittedTracks_.end(), TrackEqual(refTrack));
93  if (it==refittedTracks_.end())
94  throw cms::Exception("Vertex") << "Refitted track not found in list\n";
95  size_t pos = it - refittedTracks_.begin();
96  return tracks_[pos];
97 }
98 
99 Track Vertex::refittedTrack(const TrackBaseRef & track) const
100 {
101  if (refittedTracks_.empty())
102  throw cms::Exception("Vertex") << "No refitted tracks stored in vertex\n";
103  trackRef_iterator it = find(tracks_begin(), tracks_end(), track);
104  if (it==tracks_end()) throw cms::Exception("Vertex") << "Track not found in list\n";
105  size_t pos = it - tracks_begin();
106  return refittedTracks_[pos];
107 }
108 
110 {
111  return refittedTrack(TrackBaseRef(track));
112 }
113 
114 math::XYZTLorentzVectorD Vertex::p4(float mass,float minWeight) const
115 {
116 
118  ROOT::Math::LorentzVector<ROOT::Math::PxPyPzM4D<double> > vec;
119 
120  if(hasRefittedTracks()) {
121  for(std::vector<Track>::const_iterator iter = refittedTracks_.begin();
122  iter != refittedTracks_.end(); ++iter) {
123  if (trackWeight(originalTrack(*iter)) >=minWeight) {
124  vec.SetPx(iter->px());
125  vec.SetPy(iter->py());
126  vec.SetPz(iter->pz());
127  vec.SetM(mass);
128  sum += vec;
129  }
130  }
131  }
132  else
133  {
134  for(std::vector<reco::TrackBaseRef>::const_iterator iter = tracks_begin();
135  iter != tracks_end(); iter++) {
136  if (trackWeight(*iter) >=minWeight) {
137  vec.SetPx((*iter)->px());
138  vec.SetPy((*iter)->py());
139  vec.SetPz((*iter)->pz());
140  vec.SetM(mass);
141  sum += vec;
142  }
143  }
144  }
145  return sum;
146 }
147 
148 unsigned int Vertex::nTracks(float minWeight) const
149 {
150  int n=0;
151  if(hasRefittedTracks()) {
152  for(std::vector<Track>::const_iterator iter = refittedTracks_.begin(); iter != refittedTracks_.end(); ++iter)
153  if (trackWeight(originalTrack(*iter)) >=minWeight)
154  n++;
155  }
156  else
157  {
158  for(std::vector<reco::TrackBaseRef>::const_iterator iter = tracks_begin(); iter != tracks_end(); iter++)
159  if (trackWeight(*iter) >=minWeight)
160  n++;
161  }
162  return n;
163 }
164 
int i
Definition: DBlmapReader.cc:9
std::vector< TrackBaseRef > tracks_
reference to tracks
Definition: Vertex.h:161
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< double > > XYZTLorentzVectorD
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:15
trackRef_iterator tracks_end() const
last iterator over tracks
Definition: Vertex.cc:45
TrackBaseRef originalTrack(const Track &refTrack) const
Definition: Vertex.cc:87
Track refittedTrack(const TrackBaseRef &track) const
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:121
math::Error< dimension >::type Error
covariance error matrix (3x3)
Definition: Vertex.h:44
void fill(CovarianceMatrix &v) const
fill SMatrix
Definition: Vertex.cc:28
unsigned int index
index type
Definition: Vertex.h:50
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
edm::RefToBase< reco::Track > TrackBaseRef
persistent reference to a Track, using views
Definition: TrackFwd.h:22
float trackWeight(const TrackBaseRef &r) const
returns the weight with which a Track has contributed to the vertex-fit.
int j
Definition: DBlmapReader.cc:9
math::XYZPoint Point
point in the space
Definition: Vertex.h:40
Double32_t covariance_[size]
covariance matrix (3x3) as vector
Definition: Vertex.h:159
void add(const TrackBaseRef &r, float w=1.0)
add a reference to a Track
std::vector< float > weights_
Definition: Vertex.h:164
unsigned int nTracks(float minWeight=0.5) const
Returns the number of tracks in the vertex with weight above minWeight.
Definition: Vertex.cc:148
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:114
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector&lt;TrackRef&gt;
Definition: Vertex.h:38
void removeTracks()
Definition: Vertex.cc:66
trackRef_iterator tracks_begin() const
first iterator over tracks
Definition: Vertex.cc:40
index idx(index i, index j) const
position index
Definition: Vertex.h:170
tuple size
Write out results.
size_t tracksSize() const
number of tracks
Definition: Vertex.cc:35
bool validity_
tells wether the vertex is really valid.
Definition: Vertex.h:166