Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include "PhysicsTools/PatAlgos/interface/PATPrimaryVertexSelector.h"
0002 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0003 #include "DataFormats/TrackReco/interface/Track.h"
0004 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0005 
0006 PATPrimaryVertexSelector::PATPrimaryVertexSelector(const edm::ParameterSet& cfg, edm::ConsumesCollector&& iC)
0007     : multiplicityCut_(cfg.getParameter<unsigned int>("minMultiplicity")),
0008       ptSumCut_(cfg.getParameter<double>("minPtSum")),
0009       trackEtaCut_(cfg.getParameter<double>("maxTrackEta")),
0010       chi2Cut_(cfg.getParameter<double>("maxNormChi2")),
0011       dr2Cut_(cfg.getParameter<double>("maxDeltaR")),
0012       dzCut_(cfg.getParameter<double>("maxDeltaZ")) {
0013   // store squared cut (avoid using sqrt() for each vertex)
0014   dr2Cut_ *= dr2Cut_;
0015 }
0016 
0017 void PATPrimaryVertexSelector::select(const edm::Handle<collection>& handle,
0018                                       const edm::Event& event,
0019                                       const edm::EventSetup& setup) {
0020   selected_.clear();
0021   // FIXME: the fixed reference point should be replaced with the measured beamspot
0022   const math::XYZPoint beamSpot(0., 0., 0.);
0023   unsigned int multiplicity;
0024   double ptSum;
0025   for (collection::const_iterator iv = handle->begin(); iv != handle->end(); ++iv) {
0026     math::XYZVector displacement(iv->position() - beamSpot);
0027     if (iv->normalizedChi2() < chi2Cut_ && fabs(displacement.z()) < dzCut_ && displacement.perp2() < dr2Cut_) {
0028       getVertexVariables(*iv, multiplicity, ptSum);
0029       if (multiplicity >= multiplicityCut_ && ptSum > ptSumCut_)
0030         selected_.push_back(&*iv);
0031     }
0032   }
0033   sort(selected_.begin(), selected_.end(), *this);
0034 }
0035 
0036 bool PATPrimaryVertexSelector::operator()(const reco::Vertex* v1, const reco::Vertex* v2) const {
0037   unsigned int mult1;
0038   double ptSum1;
0039   getVertexVariables(*v1, mult1, ptSum1);
0040   unsigned int mult2;
0041   double ptSum2;
0042   getVertexVariables(*v2, mult2, ptSum2);
0043   return ptSum1 > ptSum2;
0044 }
0045 
0046 void PATPrimaryVertexSelector::getVertexVariables(const reco::Vertex& vertex,
0047                                                   unsigned int& multiplicity,
0048                                                   double& ptSum) const {
0049   multiplicity = 0;
0050   ptSum = 0.;
0051   for (reco::Vertex::trackRef_iterator it = vertex.tracks_begin(); it != vertex.tracks_end(); ++it) {
0052     if (fabs((**it).eta()) < trackEtaCut_) {
0053       ++multiplicity;
0054       ptSum += (**it).pt();
0055     }
0056   }
0057 }