CMS 3D CMS Logo

FTLRecHitProducer.cc
Go to the documentation of this file.
5 
9 
11 
13 
15 
16  public:
17  explicit FTLRecHitProducer(const edm::ParameterSet& ps);
18  ~FTLRecHitProducer() override;
19  void produce(edm::Event& evt, const edm::EventSetup& es) override;
20 
21  private:
22 
25 
26  const std::string ftlbInstance_; // instance name of barrel hits
27  const std::string ftleInstance_; // instance name of endcap hits
28 
29  std::unique_ptr<FTLRecHitAlgoBase> barrel_,endcap_;
30 };
31 
33  ftlbURecHits_( consumes<FTLUncalibratedRecHitCollection>( ps.getParameter<edm::InputTag>("barrelUncalibratedRecHits") ) ),
34  ftleURecHits_( consumes<FTLUncalibratedRecHitCollection>( ps.getParameter<edm::InputTag>("endcapUncalibratedRecHits") ) ),
35  ftlbInstance_( ps.getParameter<std::string>("BarrelHitsName") ),
36  ftleInstance_( ps.getParameter<std::string>("EndcapHitsName") ) {
37 
38  produces< FTLRecHitCollection >(ftlbInstance_);
39  produces< FTLRecHitCollection >(ftleInstance_);
40 
41  auto sumes = consumesCollector();
42 
43  const edm::ParameterSet& barrel = ps.getParameterSet("barrel");
44  const std::string& barrelAlgo = barrel.getParameter<std::string>("algoName");
45  barrel_.reset( FTLRecHitAlgoFactory::get()->create(barrelAlgo, barrel, sumes) );
46 
47  const edm::ParameterSet& endcap = ps.getParameterSet("endcap");
48  const std::string& endcapAlgo = endcap.getParameter<std::string>("algoName");
49  endcap_.reset( FTLRecHitAlgoFactory::get()->create(endcapAlgo, endcap, sumes) );
50 }
51 
53 }
54 
55 void
57 
58 
59  // tranparently get things from event setup
60  barrel_->getEventSetup(es);
61  endcap_->getEventSetup(es);
62 
63  barrel_->getEvent(evt);
64  endcap_->getEvent(evt);
65 
66  // prepare output
67  auto barrelRechits = std::make_unique<FTLRecHitCollection>();
68  auto endcapRechits = std::make_unique<FTLRecHitCollection>();
69 
71  evt.getByToken( ftlbURecHits_, hBarrel );
72  barrelRechits->reserve(hBarrel->size()/2);
73  for(const auto& uhit : *hBarrel) {
74  uint32_t flags = FTLRecHit::kGood;
75  auto rechit = std::move(barrel_->makeRecHit(uhit, flags));
76  if( flags == FTLRecHit::kGood ) barrelRechits->push_back( std::move(rechit) );
77  }
78 
80  evt.getByToken( ftleURecHits_, hEndcap );
81  endcapRechits->reserve(hEndcap->size()/2);
82  for(const auto& uhit : *hEndcap) {
83  uint32_t flags = FTLRecHit::kGood;
84  auto rechit = std::move(endcap_->makeRecHit(uhit, flags));
85  if( flags == FTLRecHit::kGood ) endcapRechits->push_back( std::move(rechit) );
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 
T getParameter(std::string const &) const
std::unique_ptr< FTLRecHitAlgoBase > endcap_
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:136
~FTLRecHitProducer() override
def create(alignables, pedeDump, additionalData, outputFile, config)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:519
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
std::vector< Variable::Flags > flags
Definition: MVATrainer.cc:135
FTLRecHitProducer(const edm::ParameterSet &ps)
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
const std::string ftleInstance_
std::unique_ptr< FTLRecHitAlgoBase > barrel_
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
const edm::EDGetTokenT< FTLUncalibratedRecHitCollection > ftlbURecHits_
ParameterSet const & getParameterSet(std::string const &) const
const std::string ftlbInstance_
HLT enums.
size_type size() const
void produce(edm::Event &evt, const edm::EventSetup &es) override
def move(src, dest)
Definition: eostools.py:510
const edm::EDGetTokenT< FTLUncalibratedRecHitCollection > ftleURecHits_
T get(const Candidate &c)
Definition: component.h:55