Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:29:02

0001 #include <memory>
0002 #include <set>
0003 
0004 #include "FWCore/Framework/interface/global/EDProducer.h"
0005 #include "FWCore/Framework/interface/Event.h"
0006 #include "FWCore/Framework/interface/MakerMacros.h"
0007 #include "FWCore/Utilities/interface/InputTag.h"
0008 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0009 
0010 #include "DataFormats/Common/interface/Handle.h"
0011 #include "DataFormats/TrackReco/interface/Track.h"
0012 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0013 #include "DataFormats/VertexReco/interface/Vertex.h"
0014 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0015 
0016 class DoubleVertexFilter : public edm::global::EDProducer<> {
0017 public:
0018   DoubleVertexFilter(const edm::ParameterSet &params);
0019 
0020   void produce(edm::StreamID, edm::Event &event, const edm::EventSetup &es) const override;
0021 
0022 private:
0023   bool trackFilter(const reco::TrackRef &track) const;
0024 
0025   edm::EDGetTokenT<reco::VertexCollection> token_primaryVertex;
0026   edm::EDGetTokenT<reco::VertexCollection> token_secondaryVertex;
0027   double maxFraction;
0028 };
0029 
0030 DoubleVertexFilter::DoubleVertexFilter(const edm::ParameterSet &params)
0031     : maxFraction(params.getParameter<double>("maxFraction")) {
0032   token_primaryVertex = consumes<reco::VertexCollection>(params.getParameter<edm::InputTag>("primaryVertices"));
0033   token_secondaryVertex = consumes<reco::VertexCollection>(params.getParameter<edm::InputTag>("secondaryVertices"));
0034   produces<reco::VertexCollection>();
0035 }
0036 
0037 static double computeSharedTracks(const reco::Vertex &pv, const reco::Vertex &sv) {
0038   std::set<reco::TrackRef> pvTracks;
0039   for (std::vector<reco::TrackBaseRef>::const_iterator iter = pv.tracks_begin(); iter != pv.tracks_end(); iter++) {
0040     if (pv.trackWeight(*iter) >= 0.5)
0041       pvTracks.insert(iter->castTo<reco::TrackRef>());
0042   }
0043 
0044   unsigned int count = 0, total = 0;
0045   for (std::vector<reco::TrackBaseRef>::const_iterator iter = sv.tracks_begin(); iter != sv.tracks_end(); iter++) {
0046     if (sv.trackWeight(*iter) >= 0.5) {
0047       total++;
0048       count += pvTracks.count(iter->castTo<reco::TrackRef>());
0049     }
0050   }
0051 
0052   return (double)count / (double)total;
0053 }
0054 
0055 void DoubleVertexFilter::produce(edm::StreamID, edm::Event &event, const edm::EventSetup &es) const {
0056   using namespace reco;
0057 
0058   edm::Handle<VertexCollection> primaryVertices;
0059   event.getByToken(token_primaryVertex, primaryVertices);
0060 
0061   edm::Handle<VertexCollection> secondaryVertices;
0062   event.getByToken(token_secondaryVertex, secondaryVertices);
0063 
0064   std::vector<reco::Vertex>::const_iterator pv = primaryVertices->begin();
0065 
0066   auto recoVertices = std::make_unique<VertexCollection>();
0067   for (std::vector<reco::Vertex>::const_iterator sv = secondaryVertices->begin(); sv != secondaryVertices->end();
0068        ++sv) {
0069     if (computeSharedTracks(*pv, *sv) > maxFraction)
0070       continue;
0071 
0072     recoVertices->push_back(*sv);
0073   }
0074 
0075   event.put(std::move(recoVertices));
0076 }
0077 
0078 DEFINE_FWK_MODULE(DoubleVertexFilter);