Go to the documentation of this file.00001 #include "DataFormats/VertexReco/interface/Vertex.h"
00002 #include <Math/GenVector/PxPyPzE4D.h>
00003 #include <Math/GenVector/PxPyPzM4D.h>
00004
00005
00006 using namespace reco;
00007 using namespace std;
00008
00009 Vertex::Vertex( const Point & p , const Error & err, double chi2, double ndof, size_t size ) :
00010 chi2_( chi2 ), ndof_( ndof ), position_( p ) {
00011 tracks_.reserve( size );
00012 index idx = 0;
00013 for( index i = 0; i < dimension; ++ i )
00014 for( index j = 0; j <= i; ++ j )
00015 covariance_[ idx ++ ] = err( i, j );
00016 validity_ = true;
00017 }
00018
00019 Vertex::Vertex( const Point & p , const Error & err) :
00020 chi2_( 0.0 ), ndof_( 0 ), position_( p ) {
00021 index idx = 0;
00022 for( index i = 0; i < dimension; ++ i )
00023 for( index j = 0; j <= i; ++ j )
00024 covariance_[ idx ++ ] = err( i, j );
00025 validity_ = true;
00026 }
00027
00028 void Vertex::fill( Error & err ) const {
00029 index idx = 0;
00030 for( index i = 0; i < dimension; ++ i )
00031 for( index j = 0; j <= i; ++ j )
00032 err( i, j ) = covariance_[ idx ++ ];
00033 }
00034
00035 size_t Vertex::tracksSize() const
00036 {
00037 return weights_.size();
00038 }
00039
00040 Vertex::trackRef_iterator Vertex::tracks_begin() const
00041 {
00042 return tracks_.begin();
00043 }
00044
00045 Vertex::trackRef_iterator Vertex::tracks_end() const
00046 {
00047
00048 return tracks_.end();
00049
00050 }
00051
00052 void Vertex::add ( const TrackBaseRef & r, float w )
00053 {
00054 tracks_.push_back ( r );
00055 weights_.push_back(w);
00056
00057 }
00058
00059 void Vertex::add ( const TrackBaseRef & r, const Track & refTrack, float w )
00060 {
00061 tracks_.push_back ( r );
00062 refittedTracks_.push_back ( refTrack );
00063 weights_.push_back(w);
00064 }
00065
00066 void Vertex::removeTracks()
00067 {
00068 weights_.clear();
00069 tracks_.clear();
00070 refittedTracks_.clear();
00071 }
00072
00073 float Vertex::trackWeight ( const TrackBaseRef & track ) const
00074 {
00075 trackRef_iterator it = find(tracks_begin(), tracks_end(), track);
00076 if (it==tracks_end()) return 0.0;
00077 size_t pos = it - tracks_begin();
00078 return weights_[pos];
00079 }
00080
00081 float Vertex::trackWeight ( const TrackRef & track ) const
00082 {
00083 return trackWeight(TrackBaseRef(track));
00084 }
00085
00086
00087 TrackBaseRef Vertex::originalTrack(const Track & refTrack) const
00088 {
00089 if (refittedTracks_.empty())
00090 throw cms::Exception("Vertex") << "No refitted tracks stored in vertex\n";
00091 std::vector<Track>::const_iterator it =
00092 find_if(refittedTracks_.begin(), refittedTracks_.end(), TrackEqual(refTrack));
00093 if (it==refittedTracks_.end())
00094 throw cms::Exception("Vertex") << "Refitted track not found in list\n";
00095 size_t pos = it - refittedTracks_.begin();
00096 return tracks_[pos];
00097 }
00098
00099 Track Vertex::refittedTrack(const TrackBaseRef & track) const
00100 {
00101 if (refittedTracks_.empty())
00102 throw cms::Exception("Vertex") << "No refitted tracks stored in vertex\n";
00103 trackRef_iterator it = find(tracks_begin(), tracks_end(), track);
00104 if (it==tracks_end()) throw cms::Exception("Vertex") << "Track not found in list\n";
00105 size_t pos = it - tracks_begin();
00106 return refittedTracks_[pos];
00107 }
00108
00109 Track Vertex::refittedTrack(const TrackRef & track) const
00110 {
00111 return refittedTrack(TrackBaseRef(track));
00112 }
00113
00114 math::XYZTLorentzVectorD Vertex::p4(float mass,float minWeight) const
00115 {
00116
00117 math::XYZTLorentzVectorD sum;
00118 ROOT::Math::LorentzVector<ROOT::Math::PxPyPzM4D<double> > vec;
00119
00120 if(hasRefittedTracks()) {
00121 for(std::vector<Track>::const_iterator iter = refittedTracks_.begin();
00122 iter != refittedTracks_.end(); ++iter) {
00123 if (trackWeight(originalTrack(*iter)) >=minWeight) {
00124 vec.SetPx(iter->px());
00125 vec.SetPy(iter->py());
00126 vec.SetPz(iter->pz());
00127 vec.SetM(mass);
00128 sum += vec;
00129 }
00130 }
00131 }
00132 else
00133 {
00134 for(std::vector<reco::TrackBaseRef>::const_iterator iter = tracks_begin();
00135 iter != tracks_end(); iter++) {
00136 if (trackWeight(*iter) >=minWeight) {
00137 vec.SetPx((*iter)->px());
00138 vec.SetPy((*iter)->py());
00139 vec.SetPz((*iter)->pz());
00140 vec.SetM(mass);
00141 sum += vec;
00142 }
00143 }
00144 }
00145 return sum;
00146 }
00147
00148 unsigned int Vertex::nTracks(float minWeight) const
00149 {
00150 int n=0;
00151 if(hasRefittedTracks()) {
00152 for(std::vector<Track>::const_iterator iter = refittedTracks_.begin(); iter != refittedTracks_.end(); ++iter)
00153 if (trackWeight(originalTrack(*iter)) >=minWeight)
00154 n++;
00155 }
00156 else
00157 {
00158 for(std::vector<reco::TrackBaseRef>::const_iterator iter = tracks_begin(); iter != tracks_end(); iter++)
00159 if (trackWeight(*iter) >=minWeight)
00160 n++;
00161 }
00162 return n;
00163 }
00164