CMS 3D CMS Logo

HGCalConcentratorProcessorSelection.cc
Go to the documentation of this file.
2 #include <limits>
3 
5 
7 
10  fixedDataSizePerHGCROC_(conf.getParameter<bool>("fixedDataSizePerHGCROC")),
11  coarsenTriggerCells_(conf.getParameter<bool>("coarsenTriggerCells")) {
13  if (selectionType == "thresholdSelect") {
15  thresholdImpl_ = std::make_unique<HGCalConcentratorThresholdImpl>(conf);
16  } else if (selectionType == "bestChoiceSelect") {
18  bestChoiceImpl_ = std::make_unique<HGCalConcentratorBestChoiceImpl>(conf);
19  } else if (selectionType == "superTriggerCellSelect") {
21  superTriggerCellImpl_ = std::make_unique<HGCalConcentratorSuperTriggerCellImpl>(conf);
22  } else if (selectionType == "mixedBestChoiceSuperTriggerCell") {
24  bestChoiceImpl_ = std::make_unique<HGCalConcentratorBestChoiceImpl>(conf);
25  superTriggerCellImpl_ = std::make_unique<HGCalConcentratorSuperTriggerCellImpl>(conf);
26  } else {
27  throw cms::Exception("HGCTriggerParameterError")
28  << "Unknown type of concentrator selection '" << selectionType << "'";
29  }
30 
32  coarsenerImpl_ = std::make_unique<HGCalConcentratorCoarsenerImpl>(conf);
33  }
34 }
35 
37  l1t::HGCalTriggerCellBxCollection& triggerCellCollOutput,
38  const edm::EventSetup& es) {
39  if (thresholdImpl_)
40  thresholdImpl_->eventSetup(es);
41  if (bestChoiceImpl_)
42  bestChoiceImpl_->eventSetup(es);
44  superTriggerCellImpl_->eventSetup(es);
45  if (coarsenerImpl_)
46  coarsenerImpl_->eventSetup(es);
48 
49  const l1t::HGCalTriggerCellBxCollection& collInput = *triggerCellCollInput;
50 
51  std::unordered_map<uint32_t, std::vector<l1t::HGCalTriggerCell>> tc_modules;
52  for (const auto& trigCell : collInput) {
53  uint32_t module = geometry_->getModuleFromTriggerCell(trigCell.detId());
54  tc_modules[module].push_back(trigCell);
55  }
56 
57  for (const auto& module_trigcell : tc_modules) {
58  std::vector<l1t::HGCalTriggerCell> trigCellVecOutput;
59  std::vector<l1t::HGCalTriggerCell> trigCellVecCoarsened;
60 
61  int thickness = triggerTools_.thicknessIndex(module_trigcell.second.at(0).detId(), true);
62 
64  coarsenerImpl_->coarsen(module_trigcell.second, trigCellVecCoarsened);
65 
66  switch (selectionType_) {
67  case thresholdSelect:
68  thresholdImpl_->select(trigCellVecCoarsened, trigCellVecOutput);
69  break;
70  case bestChoiceSelect:
71  if (triggerTools_.isEm(module_trigcell.first)) {
72  bestChoiceImpl_->select(geometry_->getLinksInModule(module_trigcell.first),
73  geometry_->getModuleSize(module_trigcell.first),
74  module_trigcell.second,
75  trigCellVecOutput);
76  } else {
77  bestChoiceImpl_->select(geometry_->getLinksInModule(module_trigcell.first),
78  geometry_->getModuleSize(module_trigcell.first),
79  trigCellVecCoarsened,
80  trigCellVecOutput);
81  }
82  break;
84  superTriggerCellImpl_->select(trigCellVecCoarsened, trigCellVecOutput);
85  break;
87  if (triggerTools_.isEm(module_trigcell.first)) {
88  bestChoiceImpl_->select(geometry_->getLinksInModule(module_trigcell.first),
89  geometry_->getModuleSize(module_trigcell.first),
90  trigCellVecCoarsened,
91  trigCellVecOutput);
92  } else {
93  superTriggerCellImpl_->select(trigCellVecCoarsened, trigCellVecOutput);
94  }
95  break;
96  default:
97  // Should not happen, selection type checked in constructor
98  break;
99  }
100 
101  } else {
102  switch (selectionType_) {
103  case thresholdSelect:
104  thresholdImpl_->select(module_trigcell.second, trigCellVecOutput);
105  break;
106  case bestChoiceSelect:
107  bestChoiceImpl_->select(geometry_->getLinksInModule(module_trigcell.first),
108  geometry_->getModuleSize(module_trigcell.first),
109  module_trigcell.second,
110  trigCellVecOutput);
111  break;
113  superTriggerCellImpl_->select(module_trigcell.second, trigCellVecOutput);
114  break;
116  if (triggerTools_.isEm(module_trigcell.first)) {
117  bestChoiceImpl_->select(geometry_->getLinksInModule(module_trigcell.first),
118  geometry_->getModuleSize(module_trigcell.first),
119  module_trigcell.second,
120  trigCellVecOutput);
121  } else {
122  superTriggerCellImpl_->select(module_trigcell.second, trigCellVecOutput);
123  }
124  break;
125  default:
126  // Should not happen, selection type checked in constructor
127  break;
128  }
129  }
130 
131  for (const auto& trigCell : trigCellVecOutput) {
132  triggerCellCollOutput.push_back(0, trigCell);
133  }
134  }
135 }
T getParameter(std::string const &) const
std::unique_ptr< HGCalConcentratorBestChoiceImpl > bestChoiceImpl_
void eventSetup(const edm::EventSetup &)
const HGCalTriggerGeometryBase * geometry_
std::unique_ptr< HGCalConcentratorSuperTriggerCellImpl > superTriggerCellImpl_
std::unique_ptr< HGCalConcentratorThresholdImpl > thresholdImpl_
HGCalConcentratorProcessorSelection(const edm::ParameterSet &conf)
virtual unsigned getLinksInModule(const unsigned module_id) const =0
void run(const edm::Handle< l1t::HGCalTriggerCellBxCollection > &triggerCellCollInput, l1t::HGCalTriggerCellBxCollection &triggerCellCollOutput, const edm::EventSetup &es) override
virtual unsigned getModuleSize(const unsigned module_id) const =0
int thicknessIndex(const DetId &, bool tc=false) const
virtual unsigned getModuleFromTriggerCell(const unsigned trigger_cell_det_id) const =0
std::unique_ptr< HGCalConcentratorCoarsenerImpl > coarsenerImpl_
bool isEm(const DetId &) const
#define DEFINE_EDM_PLUGIN(factory, type, name)
Definition: vlib.h:198
void push_back(int bx, T object)