CMS 3D CMS Logo

MuBxSelector.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 MuBxSelector(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 minNMu_;
40  std::vector<double> minMuPt_;
41  std::vector<double> maxMuEta_;
42  std::vector<int> minMuTfIndex_;
43  std::vector<int> maxMuTfIndex_;
44  std::vector<int> minMuHwQual_;
45 };
46 
48  : muonsTokenData_(consumes(iPSet.getParameter<edm::InputTag>("muonsTag"))),
49  minNMu_(iPSet.getParameter<int>("minNMu")),
50  minMuPt_(iPSet.getParameter<std::vector<double>>("minMuPt")),
51  maxMuEta_(iPSet.getParameter<std::vector<double>>("maxMuEta")),
52  minMuTfIndex_(iPSet.getParameter<std::vector<int>>("minMuTfIndex")),
53  maxMuTfIndex_(iPSet.getParameter<std::vector<int>>("maxMuTfIndex")),
54  minMuHwQual_(iPSet.getParameter<std::vector<int>>("minMuHwQual"))
55 
56 {
57  if ((minMuPt_.size() != (size_t)(size_t)minNMu_) || (maxMuEta_.size() != (size_t)minNMu_) ||
58  (minMuTfIndex_.size() != (size_t)minNMu_) || (maxMuTfIndex_.size() != (size_t)minNMu_) ||
59  (minMuHwQual_.size() != (size_t)minNMu_))
60  throw cms::Exception("MuBxSelector::MuBxSelector")
61  << "size mismatch: size of minMuPt or maxMuEta or minMuTfIndex or maxMuTfIndex or minMuHwQual != minNMu.";
62 
63  produces<std::vector<unsigned>>("SelBx").setBranchAlias("MuSelectedBx");
64 }
65 
66 // ------------ method called for each ORBIT ------------
69 
70  iEvent.getByToken(muonsTokenData_, muonsCollection);
71 
72  std::unique_ptr<std::vector<unsigned>> muBx(new std::vector<unsigned>);
73 
74  // loop over valid bunch crossings
75  for (const unsigned& bx : muonsCollection->getFilledBxs()) {
76  const auto& muons = muonsCollection->bxIterator(bx);
77 
78  // we have at least a muon
79  if (muons.size() < minNMu_)
80  continue;
81 
82  // it must be in a certain eta region with an pT and quality threshold
83  bool muCond = false;
84  int nAccMus = 0;
85  for (const auto& muon : muons) {
86  muCond = (std::abs(ugmt::fEta(muon.hwEta())) < maxMuEta_[nAccMus]) &&
87  (muon.tfMuonIndex() <= maxMuTfIndex_[nAccMus]) && (muon.tfMuonIndex() >= minMuTfIndex_[nAccMus]) &&
88  (ugmt::fPt(muon.hwPt()) >= minMuPt_[nAccMus]) && (muon.hwQual() >= minMuHwQual_[nAccMus]);
89  if (muCond)
90  nAccMus++; // found muon meeting requirements
91  if (nAccMus == minNMu_)
92  break; // found all requested muons
93  }
94 
95  if (nAccMus < minNMu_)
96  continue;
97 
98  muBx->push_back(bx);
99 
100  } // end orbit loop
101 
102  iEvent.put(std::move(muBx), "SelBx");
103 }
104 
107  desc.setUnknown();
108  descriptions.addDefault(desc);
109 }
110 
void produce(edm::Event &, const edm::EventSetup &) override
Definition: MuBxSelector.cc:67
muons
the two sets of parameters below are mutually exclusive, depending if RECO or ALCARECO is used the us...
Definition: DiMuonV_cfg.py:214
std::vector< int > minMuTfIndex_
Definition: MuBxSelector.cc:42
int iEvent
Definition: GenABIO.cc:224
void addDefault(ParameterSetDescription const &psetDescription)
MuBxSelector(const edm::ParameterSet &)
Definition: MuBxSelector.cc:47
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::vector< double > maxMuEta_
Definition: MuBxSelector.cc:41
static void fillDescriptions(edm::ConfigurationDescriptions &)
std::vector< double > minMuPt_
Definition: MuBxSelector.cc:40
edm::EDGetTokenT< OrbitCollection< l1ScoutingRun3::Muon > > muonsTokenData_
Definition: MuBxSelector.cc:36
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
std::vector< int > minMuHwQual_
Definition: MuBxSelector.cc:44
float fPt(int hwPt)
Definition: conversion.h:16
std::vector< int > maxMuTfIndex_
Definition: MuBxSelector.cc:43
HLT enums.
def move(src, dest)
Definition: eostools.py:511
float fEta(int hwEta)
Definition: conversion.h:17