CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_4_5_patch3/src/PhysicsTools/PatAlgos/src/PATPrimaryVertexSelector.cc

Go to the documentation of this file.
00001 #include "PhysicsTools/PatAlgos/interface/PATPrimaryVertexSelector.h"
00002 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00003 #include "DataFormats/TrackReco/interface/Track.h"
00004 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00005 
00006 PATPrimaryVertexSelector::PATPrimaryVertexSelector (const edm::ParameterSet& cfg) :
00007   multiplicityCut_(cfg.getParameter<unsigned int>("minMultiplicity")),
00008   ptSumCut_(cfg.getParameter<double>("minPtSum")),
00009   trackEtaCut_(cfg.getParameter<double>("maxTrackEta")),
00010   chi2Cut_(cfg.getParameter<double>("maxNormChi2")),
00011   dr2Cut_(cfg.getParameter<double>("maxDeltaR")),
00012   dzCut_(cfg.getParameter<double>("maxDeltaZ")) {
00013   // store squared cut (avoid using sqrt() for each vertex)
00014   dr2Cut_ *= dr2Cut_;
00015 }
00016 
00017 void
00018 PATPrimaryVertexSelector::select (const edm::Handle<collection>& handle, 
00019                                   const edm::Event& event,
00020                                   const edm::EventSetup& setup) {
00021   selected_.clear();
00022   // FIXME: the fixed reference point should be replaced with the measured beamspot
00023   const math::XYZPoint beamSpot(0.,0.,0.);
00024   unsigned int multiplicity;
00025   double ptSum;
00026   for ( collection::const_iterator iv=handle->begin(); iv!=handle->end(); ++iv ) {
00027     math::XYZVector displacement(iv->position()-beamSpot);
00028     if ( iv->normalizedChi2()<chi2Cut_ &&
00029          fabs(displacement.z())<dzCut_ && displacement.perp2()<dr2Cut_ ) {
00030       getVertexVariables(*iv,multiplicity,ptSum);
00031       if ( multiplicity>=multiplicityCut_ && ptSum>ptSumCut_ ) selected_.push_back(&*iv);
00032     }
00033   }
00034   sort(selected_.begin(),selected_.end(),*this);
00035 }
00036 
00037 bool
00038 PATPrimaryVertexSelector::operator() (const reco::Vertex* v1, const reco::Vertex* v2) const {
00039   unsigned int mult1;
00040   double ptSum1;
00041   getVertexVariables(*v1,mult1,ptSum1);
00042   unsigned int mult2;
00043   double ptSum2;
00044   getVertexVariables(*v2,mult2,ptSum2);
00045   return ptSum1>ptSum2;
00046 }
00047 
00048 void
00049 PATPrimaryVertexSelector::getVertexVariables (const reco::Vertex& vertex,
00050                                               unsigned int& multiplicity, double& ptSum) const {
00051   multiplicity = 0;
00052   ptSum = 0.;
00053   for(reco::Vertex::trackRef_iterator it=vertex.tracks_begin();
00054       it!=vertex.tracks_end(); ++it) {
00055     if(fabs((**it).eta())<trackEtaCut_) {
00056       ++multiplicity;
00057       ptSum += (**it).pt();
00058     }
00059   }
00060 }