CMS 3D CMS Logo

MTDRecHitProducer.cc
Go to the documentation of this file.
5 
9 
11 
19 
21 
23 public:
24  explicit MTDRecHitProducer(const edm::ParameterSet& ps);
25  ~MTDRecHitProducer() override;
26  void produce(edm::Event& evt, const edm::EventSetup& es) override;
27 
28 private:
31 
32  const std::string ftlbInstance_; // instance name of barrel hits
33  const std::string ftleInstance_; // instance name of endcap hits
34 
35  std::unique_ptr<MTDRecHitAlgoBase> barrel_, endcap_;
36 
38 };
39 
41  : ftlbURecHits_(
42  consumes<FTLUncalibratedRecHitCollection>(ps.getParameter<edm::InputTag>("barrelUncalibratedRecHits"))),
44  consumes<FTLUncalibratedRecHitCollection>(ps.getParameter<edm::InputTag>("endcapUncalibratedRecHits"))),
45  ftlbInstance_(ps.getParameter<std::string>("BarrelHitsName")),
46  ftleInstance_(ps.getParameter<std::string>("EndcapHitsName")) {
47  produces<FTLRecHitCollection>(ftlbInstance_);
48  produces<FTLRecHitCollection>(ftleInstance_);
49 
50  auto sumes = consumesCollector();
51 
52  const edm::ParameterSet& barrel = ps.getParameterSet("barrel");
53  const std::string& barrelAlgo = barrel.getParameter<std::string>("algoName");
54  barrel_ = MTDRecHitAlgoFactory::get()->create(barrelAlgo, barrel, sumes);
55 
56  const edm::ParameterSet& endcap = ps.getParameterSet("endcap");
57  const std::string& endcapAlgo = endcap.getParameter<std::string>("algoName");
58  endcap_ = MTDRecHitAlgoFactory::get()->create(endcapAlgo, endcap, sumes);
59 }
60 
62 
65  es.get<MTDDigiGeometryRecord>().get(geom);
66  geom_ = geom.product();
67 
68  // tranparently get things from event setup
69  barrel_->getEventSetup(es);
70  endcap_->getEventSetup(es);
71 
72  barrel_->getEvent(evt);
73  endcap_->getEvent(evt);
74 
75  // prepare output
76  auto barrelRechits = std::make_unique<FTLRecHitCollection>();
77  auto endcapRechits = std::make_unique<FTLRecHitCollection>();
78 
80  evt.getByToken(ftlbURecHits_, hBarrel);
81  barrelRechits->reserve(hBarrel->size() / 2);
82  for (const auto& uhit : *hBarrel) {
83  uint32_t flags = FTLRecHit::kGood;
84  auto rechit = barrel_->makeRecHit(uhit, flags);
85  if (flags == FTLRecHit::kGood)
86  barrelRechits->push_back(std::move(rechit));
87  }
88 
90  evt.getByToken(ftleURecHits_, hEndcap);
91  endcapRechits->reserve(hEndcap->size() / 2);
92  for (const auto& uhit : *hEndcap) {
93  uint32_t flags = FTLRecHit::kGood;
94  auto rechit = endcap_->makeRecHit(uhit, flags);
95  if (flags == FTLRecHit::kGood)
96  endcapRechits->push_back(std::move(rechit));
97  }
98 
99  // put the collection of recunstructed hits in the event
100  // get the orphan handles so we can make refs for the tracking rechits
101  auto barrelHandle = evt.put(std::move(barrelRechits), ftlbInstance_);
102  auto endcapHandle = evt.put(std::move(endcapRechits), ftleInstance_);
103 }
104 
T getParameter(std::string const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:131
~MTDRecHitProducer() override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
MTDRecHitProducer(const edm::ParameterSet &ps)
std::unique_ptr< MTDRecHitAlgoBase > endcap_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
const edm::EDGetTokenT< FTLUncalibratedRecHitCollection > ftleURecHits_
const MTDGeometry * geom_
const std::string ftleInstance_
const std::string ftlbInstance_
ParameterSet const & getParameterSet(std::string const &) const
void produce(edm::Event &evt, const edm::EventSetup &es) override
HLT enums.
size_type size() const
T get() const
Definition: EventSetup.h:73
std::unique_ptr< MTDRecHitAlgoBase > barrel_
const edm::EDGetTokenT< FTLUncalibratedRecHitCollection > ftlbURecHits_
T const * product() const
Definition: ESHandle.h:86
def move(src, dest)
Definition: eostools.py:511