CMS 3D CMS Logo

HGCalVFEProducer.cc
Go to the documentation of this file.
4 
10 
13 
15 
16 #include <memory>
17 
19 public:
21  ~HGCalVFEProducer() override {}
22 
23  void beginRun(const edm::Run&, const edm::EventSetup&) override;
24  void produce(edm::Event&, const edm::EventSetup&) override;
25 
26 private:
27  // inputs
31  std::unique_ptr<HGCalVFEProcessorBase> vfeProcess_;
32 };
33 
35 
37  : inputee_(consumes<HGCalDigiCollection>(conf.getParameter<edm::InputTag>("eeDigis"))),
38  inputfh_(consumes<HGCalDigiCollection>(conf.getParameter<edm::InputTag>("fhDigis"))),
39  inputbh_(consumes<HGCalDigiCollection>(conf.getParameter<edm::InputTag>("bhDigis"))),
40  triggerGeomToken_(esConsumes<HGCalTriggerGeometryBase, CaloGeometryRecord, edm::Transition::BeginRun>()) {
41  // setup VFE parameters
42  const edm::ParameterSet& vfeParamConfig = conf.getParameterSet("ProcessorParameters");
43  const std::string& vfeProcessorName = vfeParamConfig.getParameter<std::string>("ProcessorName");
44  vfeProcess_ = std::unique_ptr<HGCalVFEProcessorBase>{
45  HGCalVFEProcessorBaseFactory::get()->create(vfeProcessorName, vfeParamConfig)};
46 
47  produces<l1t::HGCalTriggerCellBxCollection>(vfeProcess_->name());
48 }
49 
50 void HGCalVFEProducer::beginRun(const edm::Run& /*run*/, const edm::EventSetup& es) {
52  vfeProcess_->setGeometry(triggerGeometry_.product());
53 }
54 
56  // Output collection
57  auto vfe_trigcell_output = std::make_unique<l1t::HGCalTriggerCellBxCollection>();
58 
59  // Input collections
63 
64  e.getByToken(inputee_, ee_digis_h);
65  e.getByToken(inputfh_, fh_digis_h);
66  e.getByToken(inputbh_, bh_digis_h);
67 
68  // Processing DigiCollections and putting the results into the HGCalTriggerCellBxCollectio
69  if (ee_digis_h.isValid()) {
70  const HGCalDigiCollection& ee_digis = *ee_digis_h;
71  vfeProcess_->run(ee_digis, *vfe_trigcell_output);
72  }
73 
74  if (fh_digis_h.isValid()) {
75  const HGCalDigiCollection& fh_digis = *fh_digis_h;
76  vfeProcess_->run(fh_digis, *vfe_trigcell_output);
77  }
78 
79  if (bh_digis_h.isValid()) {
80  const HGCalDigiCollection& bh_digis = *bh_digis_h;
81  vfeProcess_->run(bh_digis, *vfe_trigcell_output);
82  }
83 
84  // Put in the event
85  e.put(std::move(vfe_trigcell_output), vfeProcess_->name());
86 }
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
std::unique_ptr< HGCalVFEProcessorBase > vfeProcess_
ParameterSet const & getParameterSet(std::string const &) const
void produce(edm::Event &, const edm::EventSetup &) override
edm::EDGetToken inputbh_
void beginRun(const edm::Run &, const edm::EventSetup &) override
edm::ESGetToken< HGCalTriggerGeometryBase, CaloGeometryRecord > triggerGeomToken_
T const * product() const
Definition: ESHandle.h:86
HGCalVFEProducer(const edm::ParameterSet &)
Transition
Definition: Transition.h:12
edm::EDGetToken inputee_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:130
edm::EDGetToken inputfh_
bool isValid() const
Definition: HandleBase.h:70
~HGCalVFEProducer() override
HLT enums.
#define get
edm::ESHandle< HGCalTriggerGeometryBase > triggerGeometry_
def move(src, dest)
Definition: eostools.py:511
Definition: Run.h:45