File indexing completed on 2023-03-17 10:53:46
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021 #include <memory>
0022
0023
0024 #include "FWCore/Framework/interface/Frameworkfwd.h"
0025 #include "FWCore/Framework/interface/stream/EDFilter.h"
0026
0027 #include "FWCore/Framework/interface/Event.h"
0028 #include "FWCore/Framework/interface/MakerMacros.h"
0029
0030 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0031 #include "FWCore/Utilities/interface/InputTag.h"
0032 #include "DataFormats/VertexReco/interface/Vertex.h"
0033 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0034 #include "DataFormats/BTauReco/interface/SecondaryVertexTagInfo.h"
0035
0036
0037
0038
0039 class SecondaryVertexFilter : public edm::stream::EDFilter<> {
0040 public:
0041 explicit SecondaryVertexFilter(const edm::ParameterSet&);
0042 ~SecondaryVertexFilter() override;
0043
0044 private:
0045 bool filter(edm::Event&, const edm::EventSetup&) override;
0046 edm::InputTag vertexSrc;
0047 unsigned int minNumTracks;
0048 double maxAbsZ;
0049 double maxd0;
0050
0051 };
0052
0053 SecondaryVertexFilter::SecondaryVertexFilter(const edm::ParameterSet& iConfig) {
0054 vertexSrc = iConfig.getParameter<edm::InputTag>("vertexCollection");
0055 minNumTracks = iConfig.getParameter<unsigned int>("minimumNumberOfTracks");
0056 maxAbsZ = iConfig.getParameter<double>("maxAbsZ");
0057 maxd0 = iConfig.getParameter<double>("maxd0");
0058 }
0059
0060 SecondaryVertexFilter::~SecondaryVertexFilter() {}
0061
0062 bool SecondaryVertexFilter::filter(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0063 bool result = false;
0064 edm::Handle<reco::SecondaryVertexTagInfoCollection> pvHandle;
0065 iEvent.getByLabel(vertexSrc, pvHandle);
0066 const reco::SecondaryVertexTagInfoCollection& vertices = *pvHandle.product();
0067 for (reco::SecondaryVertexTagInfoCollection::const_iterator it = vertices.begin(); it != vertices.end(); ++it) {
0068 if (it->nVertices() > 0)
0069 result = true;
0070
0071
0072
0073
0074 }
0075
0076 return result;
0077 }
0078
0079
0080 DEFINE_FWK_MODULE(SecondaryVertexFilter);