Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-02-13 02:58:22

0001 #include "FWCore/Framework/interface/Event.h"
0002 #include "FWCore/Framework/interface/Frameworkfwd.h"
0003 #include "FWCore/Framework/interface/stream/EDFilter.h"
0004 #include "FWCore/Framework/interface/MakerMacros.h"
0005 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0006 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0007 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0008 #include "FWCore/Utilities/interface/EDGetToken.h"
0009 #include "FWCore/Utilities/interface/EDPutToken.h"
0010 #include "FWCore/Utilities/interface/Exception.h"
0011 #include "FWCore/Utilities/interface/InputTag.h"
0012 #include "FWCore/Utilities/interface/StreamID.h"
0013 
0014 #include <vector>
0015 #include <set>
0016 
0017 /*
0018  * Filter orbits that don't contain at least one selected BX
0019  * from a BxSelector module and produce a vector of selected BXs
0020  */
0021 class FinalBxSelector : public edm::stream::EDFilter<> {
0022 public:
0023   explicit FinalBxSelector(const edm::ParameterSet&);
0024   ~FinalBxSelector() override {}
0025   static void fillDescriptions(edm::ConfigurationDescriptions&);
0026 
0027 private:
0028   bool filter(edm::Event&, const edm::EventSetup&) override;
0029 
0030   // tokens for BX selected by each analysis
0031   std::vector<edm::EDGetTokenT<std::vector<unsigned>>> selectedBxsToken_;
0032 };
0033 
0034 FinalBxSelector::FinalBxSelector(const edm::ParameterSet& iPSet) {
0035   // get the list of selected BXs
0036   std::vector<edm::InputTag> bxLabels = iPSet.getParameter<std::vector<edm::InputTag>>("analysisLabels");
0037   for (const auto& bxLabel : bxLabels) {
0038     selectedBxsToken_.push_back(consumes<std::vector<unsigned>>(bxLabel));
0039   }
0040 
0041   produces<std::vector<unsigned>>("SelBx").setBranchAlias("SelectedBxs");
0042 }
0043 
0044 // ------------ method called for each ORBIT  ------------
0045 bool FinalBxSelector::filter(edm::Event& iEvent, const edm::EventSetup&) {
0046   bool noBxSelected = true;
0047   std::set<unsigned> uniqueBxs;
0048 
0049   for (const auto& token : selectedBxsToken_) {
0050     edm::Handle<std::vector<unsigned>> bxList;
0051     iEvent.getByToken(token, bxList);
0052 
0053     for (const unsigned& bx : *bxList) {
0054       uniqueBxs.insert(bx);
0055       noBxSelected = false;
0056     }
0057   }
0058 
0059   auto selectedBxs = std::make_unique<std::vector<unsigned>>(uniqueBxs.begin(), uniqueBxs.end());
0060   iEvent.put(std::move(selectedBxs), "SelBx");
0061 
0062   return !noBxSelected;
0063 }
0064 
0065 void FinalBxSelector::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0066   edm::ParameterSetDescription desc;
0067   desc.setUnknown();
0068   descriptions.addDefault(desc);
0069 }
0070 
0071 DEFINE_FWK_MODULE(FinalBxSelector);