CMS 3D CMS Logo

HGCalTriggerDigiProducer.cc
Go to the documentation of this file.
3 
7 
11 
13 
17 
18 #include <memory>
19 
21  public:
24 
25  void beginRun(const edm::Run&,
26  const edm::EventSetup&) override;
27  void produce(edm::Event&, const edm::EventSetup&) override;
28 
29  private:
30  // inputs
33  // algorithm containers
34  std::unique_ptr<HGCalTriggerFECodecBase> codec_;
35  std::unique_ptr<HGCalTriggerBackendProcessor> backEndProcessor_;
36 };
37 
39 
42  inputee_(consumes<HGCalDigiCollection>(conf.getParameter<edm::InputTag>("eeDigis"))),
43  inputfh_(consumes<HGCalDigiCollection>(conf.getParameter<edm::InputTag>("fhDigis"))),
44  inputbh_(consumes<HGCalDigiCollection>(conf.getParameter<edm::InputTag>("bhDigis"))),
45  backEndProcessor_(new HGCalTriggerBackendProcessor(conf.getParameterSet("BEConfiguration"),consumesCollector()) )
46 {
47  //setup FE codec
48  const edm::ParameterSet& feCodecConfig = conf.getParameterSet("FECodec");
49  const std::string& feCodecName = feCodecConfig.getParameter<std::string>("CodecName");
50  HGCalTriggerFECodecBase* codec = HGCalTriggerFECodecFactory::get()->create(feCodecName,feCodecConfig);
51  codec_.reset(codec);
52  codec_->unSetDataPayload();
53 
54  produces<l1t::HGCFETriggerDigiCollection>();
55  // register backend processor products
56  backEndProcessor_->setProduces(*this);
57 }
58 
60  const edm::EventSetup& es) {
62  codec_->setGeometry(triggerGeometry_.product());
64 
65 }
66 
68  std::unique_ptr<l1t::HGCFETriggerDigiCollection>
69  fe_output( new l1t::HGCFETriggerDigiCollection );
70 
74 
75  e.getByToken(inputee_,ee_digis_h);
76  e.getByToken(inputfh_,fh_digis_h);
77  e.getByToken(inputbh_,bh_digis_h);
78 
79  const HGCalDigiCollection& ee_digis = *ee_digis_h;
80  const HGCalDigiCollection& fh_digis = *fh_digis_h;
81  const HGCalDigiCollection& bh_digis = *bh_digis_h;
82 
83  // First find modules containing hits and prepare list of hits for each module
84  std::unordered_map<uint32_t, HGCalDigiCollection> hit_modules_ee;
85  for(const auto& eedata : ee_digis) {
86  uint32_t module = triggerGeometry_->getModuleFromCell(eedata.id());
87  if(triggerGeometry_->disconnectedModule(module)) continue;
88  auto itr_insert = hit_modules_ee.emplace(module,HGCalDigiCollection());
89  itr_insert.first->second.push_back(eedata);
90  }
91  std::unordered_map<uint32_t,HGCalDigiCollection> hit_modules_fh;
92  for(const auto& fhdata : fh_digis) {
93  uint32_t module = triggerGeometry_->getModuleFromCell(fhdata.id());
94  if(triggerGeometry_->disconnectedModule(module)) continue;
95  auto itr_insert = hit_modules_fh.emplace(module, HGCalDigiCollection());
96  itr_insert.first->second.push_back(fhdata);
97  }
98  std::unordered_map<uint32_t,HGCalDigiCollection> hit_modules_bh;
99  for(const auto& bhdata : bh_digis) {
100  if(HcalDetId(bhdata.id()).subdetId()!=HcalEndcap) continue;
101  uint32_t module = triggerGeometry_->getModuleFromCell(bhdata.id());
102  if(triggerGeometry_->disconnectedModule(module)) continue;
103  auto itr_insert = hit_modules_bh.emplace(module, HGCalDigiCollection());
104  itr_insert.first->second.push_back(bhdata);
105  }
106  // loop on modules containing hits and call front-end processing
107  // we produce one output trigger digi per module in the FE
108  fe_output->reserve(hit_modules_ee.size() + hit_modules_fh.size() + hit_modules_bh.size());
109  for( const auto& module_hits : hit_modules_ee ) {
110  fe_output->push_back(l1t::HGCFETriggerDigi());
111  l1t::HGCFETriggerDigi& digi = fe_output->back();
112  codec_->setDataPayload(module_hits.second,HGCalDigiCollection(),HGCalDigiCollection());
113  codec_->encode(digi);
114  digi.setDetId( DetId(module_hits.first) );
115  codec_->unSetDataPayload();
116  } //end loop on EE modules
117  for( const auto& module_hits : hit_modules_fh ) {
118  fe_output->push_back(l1t::HGCFETriggerDigi());
119  l1t::HGCFETriggerDigi& digi = fe_output->back();
120  codec_->setDataPayload(HGCalDigiCollection(),module_hits.second,HGCalDigiCollection());
121  codec_->encode(digi);
122  digi.setDetId( DetId(module_hits.first) );
123  codec_->unSetDataPayload();
124  } //end loop on FH modules
125  for( const auto& module_hits : hit_modules_bh ) {
126  fe_output->push_back(l1t::HGCFETriggerDigi());
127  l1t::HGCFETriggerDigi& digi = fe_output->back();
128  codec_->setDataPayload(HGCalDigiCollection(),HGCalDigiCollection(),module_hits.second);
129  codec_->encode(digi);
130  digi.setDetId( DetId(module_hits.first) );
131  codec_->unSetDataPayload();
132  } //end loop on BH modules
133 
134 
135  // get the orphan handle and fe digi collection
136  auto fe_digis_handle = e.put(std::move(fe_output));
137  auto fe_digis_coll = *fe_digis_handle;
138 
139  //now we run the emulation of the back-end processor
140  backEndProcessor_->reset();
141  backEndProcessor_->run(fe_digis_coll,es,e);
142  backEndProcessor_->putInEvent(e);
143 }
T getParameter(std::string const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:137
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
ParameterSet const & getParameterSet(ParameterSetID const &id)
std::unique_ptr< HGCalTriggerFECodecBase > codec_
void beginRun(const edm::Run &, const edm::EventSetup &) override
void setDetId(const IDTYPE &id)
HGCalTriggerDigiProducer(const edm::ParameterSet &)
edm::ESHandle< HGCalTriggerGeometryBase > triggerGeometry_
std::unique_ptr< HGCalTriggerBackendProcessor > backEndProcessor_
virtual bool disconnectedModule(const unsigned module_id) const =0
virtual unsigned getModuleFromCell(const unsigned cell_det_id) const =0
Definition: DetId.h:18
void produce(edm::Event &, const edm::EventSetup &) override
ParameterSet const & getParameterSet(std::string const &) const
HLT enums.
T get() const
Definition: EventSetup.h:63
T const * product() const
Definition: ESHandle.h:86
Definition: vlib.h:208
def move(src, dest)
Definition: eostools.py:510
T get(const Candidate &c)
Definition: component.h:55
Definition: Run.h:44
edm::SortedCollection< HGCalDataFrame > HGCalDigiCollection