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  allTrigCellsInTrigSums_(conf.getParameter<bool>("allTrigCellsInTrigSums")),
12  coarsenTriggerCells_(conf.getParameter<std::vector<unsigned>>("coarsenTriggerCells")),
13  selectionType_(kNSubDetectors_) {
14  std::vector<std::string> selectionType(conf.getParameter<std::vector<std::string>>("Method"));
16  throw cms::Exception("HGCTriggerParameterError")
17  << "Inconsistent number of sub-detectors (should be " << kNSubDetectors_ << ")";
18  }
19 
20  for (int subdet = 0; subdet < kNSubDetectors_; subdet++) {
21  if (selectionType[subdet] == "thresholdSelect") {
23  if (!thresholdImpl_)
24  thresholdImpl_ = std::make_unique<HGCalConcentratorThresholdImpl>(conf);
25  if (!trigSumImpl_)
26  trigSumImpl_ = std::make_unique<HGCalConcentratorTrigSumImpl>(conf);
27  } else if (selectionType[subdet] == "bestChoiceSelect") {
29  if (!bestChoiceImpl_)
30  bestChoiceImpl_ = std::make_unique<HGCalConcentratorBestChoiceImpl>(conf);
31  if (!trigSumImpl_)
32  trigSumImpl_ = std::make_unique<HGCalConcentratorTrigSumImpl>(conf);
33  } else if (selectionType[subdet] == "superTriggerCellSelect") {
36  superTriggerCellImpl_ = std::make_unique<HGCalConcentratorSuperTriggerCellImpl>(conf);
37  } else if (selectionType[subdet] == "autoEncoder") {
39  if (!autoEncoderImpl_)
40  autoEncoderImpl_ = std::make_unique<HGCalConcentratorAutoEncoderImpl>(conf);
41  } else if (selectionType[subdet] == "noSelection") {
42  selectionType_[subdet] = noSelection;
43  } else {
44  throw cms::Exception("HGCTriggerParameterError")
45  << "Unknown type of concentrator selection '" << selectionType[subdet] << "'";
46  }
47  }
48 
49  if (std::find(coarsenTriggerCells_.begin(), coarsenTriggerCells_.end(), true) != coarsenTriggerCells_.end() ||
51  coarsenerImpl_ = std::make_unique<HGCalConcentratorCoarsenerImpl>(conf);
52  }
53 }
54 
58  l1t::HGCalConcentratorDataBxCollection>& triggerCollOutput) {
59  if (thresholdImpl_)
60  thresholdImpl_->setGeometry(geometry());
61  if (bestChoiceImpl_)
62  bestChoiceImpl_->setGeometry(geometry());
64  superTriggerCellImpl_->setGeometry(geometry());
65  if (autoEncoderImpl_)
66  autoEncoderImpl_->setGeometry(geometry());
67  if (coarsenerImpl_)
68  coarsenerImpl_->setGeometry(geometry());
69  if (trigSumImpl_)
70  trigSumImpl_->setGeometry(geometry());
72 
73  auto& triggerCellCollOutput = std::get<0>(triggerCollOutput);
74  auto& triggerSumCollOutput = std::get<1>(triggerCollOutput);
75  auto& autoEncoderCollOutput = std::get<2>(triggerCollOutput);
76 
77  const l1t::HGCalTriggerCellBxCollection& collInput = *triggerCellCollInput;
78 
79  std::unordered_map<uint32_t, std::vector<l1t::HGCalTriggerCell>> tc_modules;
80  for (const auto& trigCell : collInput) {
81  uint32_t module = geometry()->getModuleFromTriggerCell(trigCell.detId());
82  tc_modules[module].push_back(trigCell);
83  }
84 
85  for (const auto& module_trigcell : tc_modules) {
86  std::vector<l1t::HGCalTriggerCell> trigCellVecOutput;
87  std::vector<l1t::HGCalTriggerCell> trigCellVecCoarsened;
88  std::vector<l1t::HGCalTriggerCell> trigCellVecNotSelected;
89  std::vector<l1t::HGCalTriggerSums> trigSumsVecOutput;
90  std::vector<l1t::HGCalConcentratorData> ae_EncodedLayerOutput;
91 
92  int thickness = triggerTools_.thicknessIndex(module_trigcell.second.at(0).detId());
93 
94  HGCalTriggerTools::SubDetectorType subdet = triggerTools_.getSubDetectorType(module_trigcell.second.at(0).detId());
95 
97  coarsenerImpl_->coarsen(module_trigcell.second, trigCellVecCoarsened);
98 
99  switch (selectionType_[subdet]) {
100  case thresholdSelect:
101  thresholdImpl_->select(trigCellVecCoarsened, trigCellVecOutput, trigCellVecNotSelected);
102  break;
103  case bestChoiceSelect:
104  if (triggerTools_.isEm(module_trigcell.first)) {
105  bestChoiceImpl_->select(geometry()->getLinksInModule(module_trigcell.first),
106  geometry()->getModuleSize(module_trigcell.first),
107  module_trigcell.second,
108  trigCellVecOutput,
109  trigCellVecNotSelected);
110  } else {
111  bestChoiceImpl_->select(geometry()->getLinksInModule(module_trigcell.first),
112  geometry()->getModuleSize(module_trigcell.first),
113  trigCellVecCoarsened,
114  trigCellVecOutput,
115  trigCellVecNotSelected);
116  }
117  break;
119  superTriggerCellImpl_->select(trigCellVecCoarsened, trigCellVecOutput);
120  break;
121  case autoEncoderSelect:
122  autoEncoderImpl_->select(geometry()->getLinksInModule(module_trigcell.first),
123  trigCellVecCoarsened,
124  trigCellVecOutput,
125  ae_EncodedLayerOutput);
126  break;
127  case noSelection:
128  trigCellVecOutput = trigCellVecCoarsened;
129  break;
130  default:
131  // Should not happen, selection type checked in constructor
132  break;
133  }
134 
135  } else {
136  switch (selectionType_[subdet]) {
137  case thresholdSelect:
138  thresholdImpl_->select(module_trigcell.second, trigCellVecOutput, trigCellVecNotSelected);
139  break;
140  case bestChoiceSelect:
141  bestChoiceImpl_->select(geometry()->getLinksInModule(module_trigcell.first),
142  geometry()->getModuleSize(module_trigcell.first),
143  module_trigcell.second,
144  trigCellVecOutput,
145  trigCellVecNotSelected);
146  break;
148  superTriggerCellImpl_->select(module_trigcell.second, trigCellVecOutput);
149  break;
150  case autoEncoderSelect:
151  autoEncoderImpl_->select(geometry()->getLinksInModule(module_trigcell.first),
152  module_trigcell.second,
153  trigCellVecOutput,
154  ae_EncodedLayerOutput);
155  break;
156  case noSelection:
157  trigCellVecOutput = module_trigcell.second;
158  break;
159  default:
160  // Should not happen, selection type checked in constructor
161  break;
162  }
163  }
164 
165  // trigger sum
166  if (trigSumImpl_) {
167  if (allTrigCellsInTrigSums_) { // using all TCs
168  trigSumImpl_->doSum(module_trigcell.first, module_trigcell.second, trigSumsVecOutput);
169  } else { // using only unselected TCs
170  trigSumImpl_->doSum(module_trigcell.first, trigCellVecNotSelected, trigSumsVecOutput);
171  }
172  }
173 
174  for (const auto& trigCell : trigCellVecOutput) {
175  triggerCellCollOutput.push_back(0, trigCell);
176  }
177  for (const auto& trigSums : trigSumsVecOutput) {
178  triggerSumCollOutput.push_back(0, trigSums);
179  }
180  for (const auto& aeVal : ae_EncodedLayerOutput) {
181  autoEncoderCollOutput.push_back(0, aeVal);
182  }
183  }
184 }
std::unique_ptr< HGCalConcentratorBestChoiceImpl > bestChoiceImpl_
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
int thicknessIndex(const DetId &) const
std::unique_ptr< HGCalConcentratorSuperTriggerCellImpl > superTriggerCellImpl_
std::unique_ptr< HGCalConcentratorThresholdImpl > thresholdImpl_
void run(const edm::Handle< l1t::HGCalTriggerCellBxCollection > &triggerCellCollInput, std::tuple< l1t::HGCalTriggerCellBxCollection, l1t::HGCalTriggerSumsBxCollection, l1t::HGCalConcentratorDataBxCollection > &triggerCollOutput) override
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
HGCalConcentratorProcessorSelection(const edm::ParameterSet &conf)
void setGeometry(const HGCalTriggerGeometryBase *const)
SubDetectorType getSubDetectorType(const DetId &id) const
const HGCalTriggerGeometryBase * geometry() const
virtual unsigned getModuleFromTriggerCell(const unsigned trigger_cell_det_id) const =0
bool isEm(const DetId &) const
std::unique_ptr< HGCalConcentratorAutoEncoderImpl > autoEncoderImpl_
std::unique_ptr< HGCalConcentratorCoarsenerImpl > coarsenerImpl_
#define DEFINE_EDM_PLUGIN(factory, type, name)
std::unique_ptr< HGCalConcentratorTrigSumImpl > trigSumImpl_