CMS 3D CMS Logo

JetBxSelector.cc
Go to the documentation of this file.
14 
15 // L1 scouting
19 
20 #include <memory>
21 #include <utility>
22 #include <vector>
23 
24 using namespace l1ScoutingRun3;
25 
27 public:
28  explicit JetBxSelector(const edm::ParameterSet&);
31 
32 private:
33  void produce(edm::Event&, const edm::EventSetup&) override;
34 
35  // tokens for scouting data
37 
38  // SELECTION THRESHOLDS
39  int minNJet_;
40  std::vector<double> minJetEt_;
41  std::vector<double> maxJetEta_;
42 };
43 
45  : jetsTokenData_(consumes(iPSet.getParameter<edm::InputTag>("jetsTag"))),
46  minNJet_(iPSet.getParameter<int>("minNJet")),
47  minJetEt_(iPSet.getParameter<std::vector<double>>("minJetEt")),
48  maxJetEta_(iPSet.getParameter<std::vector<double>>("maxJetEta"))
49 
50 {
51  if ((minJetEt_.size() != (size_t)minNJet_) || (maxJetEta_.size() != (size_t)minNJet_))
52  throw cms::Exception("JetBxSelector::JetBxSelector") << "size mismatch: size of minJetEt or maxJetEta != minNJet.";
53 
54  produces<std::vector<unsigned>>("SelBx").setBranchAlias("JetSelectedBx");
55 }
56 
57 // ------------ method called for each ORBIT ------------
60 
61  iEvent.getByToken(jetsTokenData_, jetsCollection);
62 
63  std::unique_ptr<std::vector<unsigned>> jetBx(new std::vector<unsigned>);
64 
65  // loop over valid bunch crossings
66  for (const unsigned& bx : jetsCollection->getFilledBxs()) {
67  const auto& jets = jetsCollection->bxIterator(bx);
68 
69  // we have at least N jets
70  if (jets.size() < minNJet_)
71  continue;
72 
73  // it must be in a certain eta region with an pT and quality threshold
74  bool jetCond = false;
75  int nAccJets = 0;
76  for (const auto& jet : jets) {
77  jetCond = (std::abs(demux::fEta(jet.hwEta())) < maxJetEta_[nAccJets]) &&
78  (demux::fEt(jet.hwEt()) >= minJetEt_[nAccJets]);
79  if (jetCond)
80  nAccJets++; // found jet meeting requirements
81  if (nAccJets == minNJet_)
82  break; // found all requested jets
83  }
84 
85  if (nAccJets < minNJet_)
86  continue;
87 
88  jetBx->push_back(bx);
89 
90  } // end orbit loop
91 
92  iEvent.put(std::move(jetBx), "SelBx");
93 }
94 
97  desc.setUnknown();
98  descriptions.addDefault(desc);
99 }
100 
float fEt(int hwEt)
Definition: conversion.h:29
std::vector< double > maxJetEta_
int iEvent
Definition: GenABIO.cc:224
void addDefault(ParameterSetDescription const &psetDescription)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
static void fillDescriptions(edm::ConfigurationDescriptions &)
void produce(edm::Event &, const edm::EventSetup &) override
JetBxSelector(const edm::ParameterSet &)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
HLT enums.
def move(src, dest)
Definition: eostools.py:511
float fEta(int hwEta)
Definition: conversion.h:17
std::vector< double > minJetEt_
edm::EDGetTokenT< OrbitCollection< l1ScoutingRun3::Jet > > jetsTokenData_