File indexing completed on 2024-04-06 12:24:33
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020 #include <memory>
0021
0022
0023 #include "FWCore/Framework/interface/Frameworkfwd.h"
0024 #include "FWCore/Framework/interface/stream/EDFilter.h"
0025
0026 #include "FWCore/Framework/interface/Event.h"
0027 #include "FWCore/Framework/interface/MakerMacros.h"
0028
0029 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0030 #include "FWCore/Utilities/interface/InputTag.h"
0031 #include "DataFormats/Candidate/interface/VertexCompositePtrCandidate.h"
0032 #include "DataFormats/VertexReco/interface/Vertex.h"
0033 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0034
0035 #include "RecoBTag/SecondaryVertex/interface/TemplatedSecondaryVertex.h"
0036 #include "RecoBTag/SecondaryVertex/interface/VertexFilter.h"
0037
0038
0039
0040
0041 template <typename VTX>
0042 class BVertexFilterT : public edm::stream::EDFilter<> {
0043 public:
0044 explicit BVertexFilterT(const edm::ParameterSet&);
0045 ~BVertexFilterT() override;
0046
0047 private:
0048 bool filter(edm::Event&, const edm::EventSetup&) override;
0049 edm::EDGetTokenT<reco::VertexCollection> token_primaryVertex;
0050 edm::EDGetTokenT<edm::View<VTX>> token_secondaryVertex;
0051 reco::VertexFilter svFilter;
0052 bool useVertexKinematicAsJetAxis;
0053 int minVertices;
0054 };
0055
0056 template <typename VTX>
0057 BVertexFilterT<VTX>::BVertexFilterT(const edm::ParameterSet& params)
0058 : svFilter(params.getParameter<edm::ParameterSet>("vertexFilter")),
0059 useVertexKinematicAsJetAxis(params.getParameter<bool>("useVertexKinematicAsJetAxis")),
0060 minVertices(params.getParameter<int>("minVertices"))
0061
0062 {
0063 token_primaryVertex = consumes<reco::VertexCollection>(params.getParameter<edm::InputTag>("primaryVertices"));
0064 token_secondaryVertex = consumes<edm::View<VTX>>(params.getParameter<edm::InputTag>("secondaryVertices"));
0065 produces<std::vector<VTX>>();
0066 }
0067
0068 template <typename VTX>
0069 BVertexFilterT<VTX>::~BVertexFilterT() {}
0070
0071 template <typename VTX>
0072 bool BVertexFilterT<VTX>::filter(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0073 int count = 0;
0074 edm::Handle<reco::VertexCollection> pvHandle;
0075 iEvent.getByToken(token_primaryVertex, pvHandle);
0076 edm::Handle<edm::View<VTX>> svHandle;
0077 iEvent.getByToken(token_secondaryVertex, svHandle);
0078
0079 auto recoVertices = std::make_unique<std::vector<VTX>>();
0080
0081 if (!pvHandle->empty()) {
0082 const reco::Vertex& primary = (*pvHandle.product())[0];
0083 const edm::View<VTX>& vertices = *svHandle.product();
0084
0085 if (!primary.isFake()) {
0086 for (typename edm::View<VTX>::const_iterator it = vertices.begin(); it != vertices.end(); ++it) {
0087 GlobalVector axis(0, 0, 0);
0088 if (useVertexKinematicAsJetAxis)
0089 axis = GlobalVector(it->p4().X(), it->p4().Y(), it->p4().Z());
0090 if (svFilter(primary, reco::TemplatedSecondaryVertex<VTX>(primary, *it, axis, true), axis)) {
0091 count++;
0092 recoVertices->push_back(*it);
0093 }
0094 }
0095 }
0096 }
0097 iEvent.put(std::move(recoVertices));
0098
0099 return (count >= minVertices);
0100 }
0101
0102
0103 typedef BVertexFilterT<reco::Vertex> BVertexFilter;
0104 typedef BVertexFilterT<reco::VertexCompositePtrCandidate> CandidateBVertexFilter;
0105
0106
0107 DEFINE_FWK_MODULE(BVertexFilter);
0108 DEFINE_FWK_MODULE(CandidateBVertexFilter);