CMS 3D CMS Logo

MTDUncalibratedRecHitProducer.cc
Go to the documentation of this file.
5 
9 
11 
13 
15 public:
18  void produce(edm::Event& evt, const edm::EventSetup& es) override;
19 
20 private:
21  const edm::EDGetTokenT<BTLDigiCollection> ftlbDigis_; // collection of BTL digis
22  const edm::EDGetTokenT<ETLDigiCollection> ftleDigis_; // collection of ETL digis
23 
24  const std::string ftlbInstance_; // instance name of barrel hits
25  const std::string ftleInstance_; // instance name of endcap hits
26 
27  std::unique_ptr<BTLUncalibratedRecHitAlgoBase> barrel_;
28  std::unique_ptr<ETLUncalibratedRecHitAlgoBase> endcap_;
29 };
30 
32  : ftlbDigis_(consumes<BTLDigiCollection>(ps.getParameter<edm::InputTag>("barrelDigis"))),
33  ftleDigis_(consumes<ETLDigiCollection>(ps.getParameter<edm::InputTag>("endcapDigis"))),
34  ftlbInstance_(ps.getParameter<std::string>("BarrelHitsName")),
35  ftleInstance_(ps.getParameter<std::string>("EndcapHitsName")) {
36  produces<FTLUncalibratedRecHitCollection>(ftlbInstance_);
37  produces<FTLUncalibratedRecHitCollection>(ftleInstance_);
38 
39  auto sumes = consumesCollector();
40 
41  const edm::ParameterSet& barrel = ps.getParameterSet("barrel");
42  const std::string& barrelAlgo = barrel.getParameter<std::string>("algoName");
43  barrel_ = std::unique_ptr<BTLUncalibratedRecHitAlgoBase>{
44  BTLUncalibratedRecHitAlgoFactory::get()->create(barrelAlgo, barrel, sumes)};
45 
46  const edm::ParameterSet& endcap = ps.getParameterSet("endcap");
47  const std::string& endcapAlgo = endcap.getParameter<std::string>("algoName");
48  endcap_ = std::unique_ptr<ETLUncalibratedRecHitAlgoBase>{
49  ETLUncalibratedRecHitAlgoFactory::get()->create(endcapAlgo, endcap, sumes)};
50 }
51 
53 
55  // tranparently get things from event setup
56  barrel_->getEventSetup(es);
57  endcap_->getEventSetup(es);
58 
59  barrel_->getEvent(evt);
60  endcap_->getEvent(evt);
61 
62  // prepare output
63  auto barrelRechits = std::make_unique<FTLUncalibratedRecHitCollection>();
64  auto endcapRechits = std::make_unique<FTLUncalibratedRecHitCollection>();
65 
67  evt.getByToken(ftlbDigis_, hBarrel);
68  if (hBarrel.isValid()) {
69  barrelRechits->reserve(hBarrel->size() / 2);
70  for (const auto& digi : *hBarrel) {
71  barrelRechits->emplace_back(barrel_->makeRecHit(digi));
72  }
73  } else {
74  edm::LogWarning("MTDReco") << "MTDUncalibratedRecHitProducer: Missing BTL Digi Collection";
75  }
76 
78  evt.getByToken(ftleDigis_, hEndcap);
79  if (hEndcap.isValid()) {
80  endcapRechits->reserve(hEndcap->size() / 2);
81  for (const auto& digi : *hEndcap) {
82  endcapRechits->emplace_back(endcap_->makeRecHit(digi));
83  }
84  } else {
85  edm::LogWarning("MTDReco") << "MTDUncalibratedRecHitProducer: Missing ETL Digi Collection";
86  }
87 
88  // put the collection of recunstructed hits in the event
89  evt.put(std::move(barrelRechits), ftlbInstance_);
90  evt.put(std::move(endcapRechits), ftleInstance_);
91 }
92 
const edm::EDGetTokenT< ETLDigiCollection > ftleDigis_
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
size_type size() const
ParameterSet const & getParameterSet(std::string const &) const
MTDUncalibratedRecHitProducer(const edm::ParameterSet &ps)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:536
std::unique_ptr< ETLUncalibratedRecHitAlgoBase > endcap_
const edm::EDGetTokenT< BTLDigiCollection > ftlbDigis_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void produce(edm::Event &evt, const edm::EventSetup &es) override
bool isValid() const
Definition: HandleBase.h:70
HLT enums.
std::unique_ptr< BTLUncalibratedRecHitAlgoBase > barrel_
#define get
Log< level::Warning, false > LogWarning
def move(src, dest)
Definition: eostools.py:511