File indexing completed on 2025-01-08 03:36:33
0001 #include <memory>
0002 #include <set>
0003
0004 #include "DataFormats/Candidate/interface/VertexCompositePtrCandidate.h"
0005 #include "DataFormats/Common/interface/Handle.h"
0006 #include "DataFormats/TrackReco/interface/Track.h"
0007 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0008 #include "DataFormats/VertexReco/interface/Vertex.h"
0009 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0010 #include "FWCore/Framework/interface/Event.h"
0011 #include "FWCore/Framework/interface/MakerMacros.h"
0012 #include "FWCore/Framework/interface/stream/EDProducer.h"
0013 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0014 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0015 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0016 #include "FWCore/Utilities/interface/InputTag.h"
0017 #include "RecoVertex/VertexPrimitives/interface/ConvertToFromReco.h"
0018 #include "RecoVertex/VertexPrimitives/interface/VertexState.h"
0019 #include "RecoVertex/VertexTools/interface/SharedTracks.h"
0020 #include "RecoVertex/VertexTools/interface/VertexDistance3D.h"
0021
0022 template <class VTX>
0023 class TemplatedVertexMerger : public edm::stream::EDProducer<> {
0024 public:
0025 typedef std::vector<VTX> Product;
0026 TemplatedVertexMerger(const edm::ParameterSet ¶ms);
0027
0028 static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
0029 void produce(edm::Event &event, const edm::EventSetup &es) override;
0030
0031 private:
0032 bool trackFilter(const reco::TrackRef &track) const;
0033
0034 edm::EDGetTokenT<Product> token_secondaryVertex;
0035 double maxFraction;
0036 double minSignificance;
0037 };
0038
0039 template <class VTX>
0040 TemplatedVertexMerger<VTX>::TemplatedVertexMerger(const edm::ParameterSet ¶ms)
0041 : maxFraction(params.getParameter<double>("maxFraction")),
0042 minSignificance(params.getParameter<double>("minSignificance")) {
0043 token_secondaryVertex = consumes<Product>(params.getParameter<edm::InputTag>("secondaryVertices"));
0044 produces<Product>();
0045 }
0046
0047 template <class VTX>
0048 void TemplatedVertexMerger<VTX>::produce(edm::Event &event, const edm::EventSetup &es) {
0049 using namespace reco;
0050
0051 edm::Handle<Product> secondaryVertices;
0052 event.getByToken(token_secondaryVertex, secondaryVertices);
0053
0054 VertexDistance3D dist;
0055 auto recoVertices = std::make_unique<Product>();
0056 for (typename Product::const_iterator sv = secondaryVertices->begin(); sv != secondaryVertices->end(); ++sv) {
0057 recoVertices->push_back(*sv);
0058 }
0059 for (typename Product::iterator sv = recoVertices->begin(); sv != recoVertices->end(); ++sv) {
0060 bool shared = false;
0061 VertexState s1(RecoVertex::convertPos(sv->position()), RecoVertex::convertError(sv->error()));
0062 for (typename Product::iterator sv2 = recoVertices->begin(); sv2 != recoVertices->end(); ++sv2) {
0063 VertexState s2(RecoVertex::convertPos(sv2->position()), RecoVertex::convertError(sv2->error()));
0064 double fr = vertexTools::computeSharedTracks(*sv2, *sv);
0065
0066
0067 if (fr > maxFraction && dist.distance(s1, s2).significance() < minSignificance && sv - sv2 != 0 &&
0068 fr >= vertexTools::computeSharedTracks(*sv, *sv2)) {
0069 shared = true;
0070
0071 }
0072 }
0073 if (shared) {
0074 sv = recoVertices->erase(sv) - 1;
0075 }
0076
0077 }
0078
0079 event.put(std::move(recoVertices));
0080 }
0081
0082 template <class VTX>
0083 void TemplatedVertexMerger<VTX>::fillDescriptions(edm::ConfigurationDescriptions &descriptions) {
0084 edm::ParameterSetDescription desc;
0085 desc.add<double>("maxFraction", 0.7);
0086 desc.add<double>("minSignificance", 2);
0087 desc.add<edm::InputTag>("secondaryVertices", edm::InputTag("inclusiveVertexFinder"));
0088 descriptions.addWithDefaultLabel(desc);
0089 }
0090
0091 typedef TemplatedVertexMerger<reco::Vertex> VertexMerger;
0092 typedef TemplatedVertexMerger<reco::VertexCompositePtrCandidate> CandidateVertexMerger;
0093
0094 DEFINE_FWK_MODULE(VertexMerger);
0095 DEFINE_FWK_MODULE(CandidateVertexMerger);