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/EDProducer.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 // L1 scouting
0015 #include "DataFormats/L1Scouting/interface/L1ScoutingCalo.h"
0016 #include "DataFormats/L1Scouting/interface/OrbitCollection.h"
0017 #include "L1TriggerScouting/Utilities/interface/conversion.h"
0018 
0019 #include <memory>
0020 #include <utility>
0021 #include <vector>
0022 
0023 using namespace l1ScoutingRun3;
0024 
0025 class JetBxSelector : public edm::stream::EDProducer<> {
0026 public:
0027   explicit JetBxSelector(const edm::ParameterSet&);
0028   ~JetBxSelector() override {}
0029   static void fillDescriptions(edm::ConfigurationDescriptions&);
0030 
0031 private:
0032   void produce(edm::Event&, const edm::EventSetup&) override;
0033 
0034   // tokens for scouting data
0035   edm::EDGetTokenT<OrbitCollection<l1ScoutingRun3::Jet>> jetsTokenData_;
0036 
0037   // SELECTION THRESHOLDS
0038   int minNJet_;
0039   std::vector<double> minJetEt_;
0040   std::vector<double> maxJetEta_;
0041 };
0042 
0043 JetBxSelector::JetBxSelector(const edm::ParameterSet& iPSet)
0044     : jetsTokenData_(consumes(iPSet.getParameter<edm::InputTag>("jetsTag"))),
0045       minNJet_(iPSet.getParameter<int>("minNJet")),
0046       minJetEt_(iPSet.getParameter<std::vector<double>>("minJetEt")),
0047       maxJetEta_(iPSet.getParameter<std::vector<double>>("maxJetEta"))
0048 
0049 {
0050   if ((minJetEt_.size() != (size_t)minNJet_) || (maxJetEta_.size() != (size_t)minNJet_))
0051     throw cms::Exception("JetBxSelector::JetBxSelector") << "size mismatch: size of minJetEt or maxJetEta != minNJet.";
0052 
0053   produces<std::vector<unsigned>>("SelBx").setBranchAlias("JetSelectedBx");
0054 }
0055 
0056 // ------------ method called for each ORBIT  ------------
0057 void JetBxSelector::produce(edm::Event& iEvent, const edm::EventSetup&) {
0058   edm::Handle<OrbitCollection<l1ScoutingRun3::Jet>> jetsCollection;
0059 
0060   iEvent.getByToken(jetsTokenData_, jetsCollection);
0061 
0062   std::unique_ptr<std::vector<unsigned>> jetBx(new std::vector<unsigned>);
0063 
0064   // loop over valid bunch crossings
0065   for (const unsigned& bx : jetsCollection->getFilledBxs()) {
0066     const auto& jets = jetsCollection->bxIterator(bx);
0067 
0068     // we have at least N jets
0069     if (std::ssize(jets) < minNJet_)
0070       continue;
0071 
0072     // it must be in a certain eta region with an pT and quality threshold
0073     bool jetCond = false;
0074     int nAccJets = 0;
0075     for (const auto& jet : jets) {
0076       jetCond = (std::abs(demux::fEta(jet.hwEta())) < maxJetEta_[nAccJets]) &&
0077                 (demux::fEt(jet.hwEt()) >= minJetEt_[nAccJets]);
0078       if (jetCond)
0079         nAccJets++;  // found jet meeting requirements
0080       if (nAccJets == minNJet_)
0081         break;  // found all requested jets
0082     }
0083 
0084     if (nAccJets < minNJet_)
0085       continue;
0086 
0087     jetBx->push_back(bx);
0088 
0089   }  // end orbit loop
0090 
0091   iEvent.put(std::move(jetBx), "SelBx");
0092 }
0093 
0094 void JetBxSelector::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0095   edm::ParameterSetDescription desc;
0096   desc.setUnknown();
0097   descriptions.addDefault(desc);
0098 }
0099 
0100 DEFINE_FWK_MODULE(JetBxSelector);