Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-12-20 03:13:59

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