Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:05:23

0001 #include "DataFormats/VertexReco/interface/Vertex.h"
0002 #include <Math/GenVector/PxPyPzE4D.h>
0003 #include <Math/GenVector/PxPyPzM4D.h>
0004 
0005 using namespace reco;
0006 using namespace std;
0007 
0008 Vertex::Vertex(const Point& p, const Error& err, double chi2, double ndof, size_t size)
0009     : chi2_(chi2), ndof_(ndof), position_(p), time_(0.) {
0010   tracks_.reserve(size);
0011   index idx = 0;
0012   for (index i = 0; i < dimension4D; ++i) {
0013     for (index j = 0; j <= i; ++j) {
0014       if (i == dimension || j == dimension) {
0015         covariance_[idx++] = 0.0;
0016       } else {
0017         covariance_[idx++] = err(i, j);
0018       }
0019     }
0020   }
0021   validity_ = true;
0022 }
0023 
0024 Vertex::Vertex(const Point& p, const Error4D& err, double time, double chi2, double ndof, size_t size)
0025     : chi2_(chi2), ndof_(ndof), position_(p), time_(time) {
0026   tracks_.reserve(size4D);
0027   index idx = 0;
0028   for (index i = 0; i < dimension4D; ++i)
0029     for (index j = 0; j <= i; ++j)
0030       covariance_[idx++] = err(i, j);
0031   validity_ = true;
0032 }
0033 
0034 Vertex::Vertex(const Point& p, const Error& err) : chi2_(0.0), ndof_(0), position_(p), time_(0.) {
0035   index idx = 0;
0036   for (index i = 0; i < dimension4D; ++i) {
0037     for (index j = 0; j <= i; ++j) {
0038       if (i == dimension || j == dimension) {
0039         covariance_[idx++] = 0.0;
0040       } else {
0041         covariance_[idx++] = err(i, j);
0042       }
0043     }
0044   }
0045   validity_ = true;
0046 }
0047 
0048 Vertex::Vertex(const Point& p, const Error4D& err, double time) : chi2_(0.0), ndof_(0), position_(p), time_(time) {
0049   index idx = 0;
0050   for (index i = 0; i < dimension + 1; ++i)
0051     for (index j = 0; j <= i; ++j)
0052       covariance_[idx++] = err(i, j);
0053   validity_ = true;
0054 }
0055 
0056 void Vertex::fill(Error& err) const {
0057   Error4D temp;
0058   fill(temp);
0059   err = temp.Sub<Error>(0, 0);
0060 }
0061 
0062 void Vertex::fill(Error4D& err) const {
0063   index idx = 0;
0064   for (index i = 0; i < dimension4D; ++i)
0065     for (index j = 0; j <= i; ++j)
0066       err(i, j) = covariance_[idx++];
0067 }
0068 
0069 void Vertex::add(const TrackBaseRef& r, const Track& refTrack, float w) {
0070   tracks_.push_back(r);
0071   refittedTracks_.push_back(refTrack);
0072   weights_.push_back(w * 255);
0073 }
0074 
0075 void Vertex::removeTracks() {
0076   weights_.clear();
0077   tracks_.clear();
0078   refittedTracks_.clear();
0079 }
0080 
0081 TrackBaseRef Vertex::originalTrack(const Track& refTrack) const {
0082   if (refittedTracks_.empty())
0083     throw cms::Exception("Vertex") << "No refitted tracks stored in vertex\n";
0084   std::vector<Track>::const_iterator it = find_if(refittedTracks_.begin(), refittedTracks_.end(), TrackEqual(refTrack));
0085   if (it == refittedTracks_.end())
0086     throw cms::Exception("Vertex") << "Refitted track not found in list.\n pt used for comparison: " << refTrack.pt();
0087   size_t pos = it - refittedTracks_.begin();
0088   return tracks_[pos];
0089 }
0090 
0091 Track Vertex::refittedTrack(const TrackBaseRef& track) const {
0092   if (refittedTracks_.empty())
0093     throw cms::Exception("Vertex") << "No refitted tracks stored in vertex\n";
0094   trackRef_iterator it = find(tracks_begin(), tracks_end(), track);
0095   if (it == tracks_end())
0096     throw cms::Exception("Vertex") << "Track not found in list\n";
0097   size_t pos = it - tracks_begin();
0098   return refittedTracks_[pos];
0099 }
0100 
0101 Track Vertex::refittedTrack(const TrackRef& track) const { return refittedTrack(TrackBaseRef(track)); }
0102 
0103 math::XYZTLorentzVectorD Vertex::p4(float mass, float minWeight) const {
0104   math::XYZTLorentzVectorD sum;
0105   ROOT::Math::LorentzVector<ROOT::Math::PxPyPzM4D<double> > vec;
0106 
0107   if (hasRefittedTracks()) {
0108     for (std::vector<Track>::const_iterator iter = refittedTracks_.begin(); iter != refittedTracks_.end(); ++iter) {
0109       if (trackWeight(originalTrack(*iter)) >= minWeight) {
0110         vec.SetPx(iter->px());
0111         vec.SetPy(iter->py());
0112         vec.SetPz(iter->pz());
0113         vec.SetM(mass);
0114         sum += vec;
0115       }
0116     }
0117   } else {
0118     for (std::vector<reco::TrackBaseRef>::const_iterator iter = tracks_begin(); iter != tracks_end(); iter++) {
0119       if (trackWeight(*iter) >= minWeight) {
0120         vec.SetPx((*iter)->px());
0121         vec.SetPy((*iter)->py());
0122         vec.SetPz((*iter)->pz());
0123         vec.SetM(mass);
0124         sum += vec;
0125       }
0126     }
0127   }
0128   return sum;
0129 }
0130 
0131 unsigned int Vertex::nTracks(float minWeight) const {
0132   int n = 0;
0133   if (hasRefittedTracks()) {
0134     for (std::vector<Track>::const_iterator iter = refittedTracks_.begin(); iter != refittedTracks_.end(); ++iter)
0135       if (trackWeight(originalTrack(*iter)) >= minWeight)
0136         n++;
0137   } else {
0138     for (std::vector<reco::TrackBaseRef>::const_iterator iter = tracks_begin(); iter != tracks_end(); iter++)
0139       if (trackWeight(*iter) >= minWeight)
0140         n++;
0141   }
0142   return n;
0143 }