Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-10-08 05:11:31

0001 #include <memory>
0002 #include <vector>
0003 #include "DataFormats/Common/interface/Handle.h"
0004 #include "DataFormats/Common/interface/View.h"
0005 #include "DataFormats/Common/interface/Ref.h"
0006 #include "DataFormats/Candidate/interface/VertexCompositeCandidate.h"
0007 #include "FWCore/Framework/interface/stream/EDFilter.h"
0008 #include "FWCore/Framework/interface/Event.h"
0009 #include "FWCore/Framework/interface/MakerMacros.h"
0010 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0011 #include "FWCore/Utilities/interface/InputTag.h"
0012 
0013 class V0EventSelector : public edm::stream::EDFilter<> {
0014 public:
0015   explicit V0EventSelector(const edm::ParameterSet&);
0016   ~V0EventSelector() override = default;
0017 
0018   bool filter(edm::Event&, const edm::EventSetup&) override;
0019   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0020 
0021 private:
0022   const edm::EDGetTokenT<reco::VertexCompositeCandidateCollection> vccToken_;
0023   const unsigned int minNumCandidates_;
0024   const double massMin_;
0025   const double massMax_;
0026 };
0027 
0028 V0EventSelector::V0EventSelector(const edm::ParameterSet& iConfig)
0029     : vccToken_{consumes<reco::VertexCompositeCandidateCollection>(
0030           iConfig.getParameter<edm::InputTag>("vertexCompositeCandidates"))},
0031       minNumCandidates_{iConfig.getParameter<unsigned int>("minCandidates")},
0032       massMin_{iConfig.getParameter<double>("massMin")},
0033       massMax_{iConfig.getParameter<double>("massMax")} {
0034   produces<reco::VertexCompositeCandidateCollection>();
0035 }
0036 
0037 bool V0EventSelector::filter(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0038   edm::Handle<reco::VertexCompositeCandidateCollection> vccHandle;
0039   iEvent.getByToken(vccToken_, vccHandle);
0040   auto filteredVCC = std::make_unique<reco::VertexCompositeCandidateCollection>();
0041 
0042   // early return if the input collection is empty
0043   if (!vccHandle.isValid()) {
0044     iEvent.put(std::move(filteredVCC));
0045     return false;
0046   }
0047 
0048   for (const auto& vcc : *vccHandle) {
0049     if (vcc.mass() >= massMin_ && vcc.mass() <= massMax_) {
0050       filteredVCC->push_back(vcc);
0051     }
0052   }
0053 
0054   bool pass = filteredVCC->size() >= minNumCandidates_;
0055   iEvent.put(std::move(filteredVCC));
0056 
0057   return pass;
0058 }
0059 
0060 void V0EventSelector::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0061   edm::ParameterSetDescription desc;
0062   desc.add<edm::InputTag>("vertexCompositeCandidates", edm::InputTag("generalV0Candidates:Kshort"));
0063   desc.add<unsigned int>("minCandidates", 1)->setComment("Minimum number of candidates required");
0064   desc.add<double>("massMin", 0.0)->setComment("Minimum mass in GeV");
0065   desc.add<double>("massMax", std::numeric_limits<double>::max())->setComment("Maximum mass in GeV");
0066   descriptions.addWithDefaultLabel(desc);
0067 }
0068 
0069 // Define this module as a plug-in
0070 DEFINE_FWK_MODULE(V0EventSelector);