CMS 3D CMS Logo

RecHitMapProducer.cc
Go to the documentation of this file.
1 // user include files
2 #include <unordered_map>
3 
7 
13 
16 
18 public:
20  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
21 
22  void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
23 
24 private:
31  bool hgcalOnly_;
32 };
33 
35 
36 using DetIdRecHitMap = std::unordered_map<DetId, const unsigned int>;
37 
39  : hits_ee_token_(consumes<HGCRecHitCollection>(ps.getParameter<edm::InputTag>("EEInput"))),
40  hits_fh_token_(consumes<HGCRecHitCollection>(ps.getParameter<edm::InputTag>("FHInput"))),
41  hits_bh_token_(consumes<HGCRecHitCollection>(ps.getParameter<edm::InputTag>("BHInput"))),
42  hits_eb_token_(consumes<reco::PFRecHitCollection>(ps.getParameter<edm::InputTag>("EBInput"))),
43  hits_hb_token_(consumes<reco::PFRecHitCollection>(ps.getParameter<edm::InputTag>("HBInput"))),
44  hits_ho_token_(consumes<reco::PFRecHitCollection>(ps.getParameter<edm::InputTag>("HOInput"))),
45  hgcalOnly_(ps.getParameter<bool>("hgcalOnly")) {
46  produces<DetIdRecHitMap>("hgcalRecHitMap");
47  if (!hgcalOnly_)
48  produces<DetIdRecHitMap>("barrelRecHitMap");
49 }
50 
53  desc.add<edm::InputTag>("EEInput", {"HGCalRecHit", "HGCEERecHits"});
54  desc.add<edm::InputTag>("FHInput", {"HGCalRecHit", "HGCHEFRecHits"});
55  desc.add<edm::InputTag>("BHInput", {"HGCalRecHit", "HGCHEBRecHits"});
56  desc.add<edm::InputTag>("EBInput", {"particleFlowRecHitECAL", ""});
57  desc.add<edm::InputTag>("HBInput", {"particleFlowRecHitHBHE", ""});
58  desc.add<edm::InputTag>("HOInput", {"particleFlowRecHitHO", ""});
59  desc.add<bool>("hgcalOnly", true);
60  descriptions.add("recHitMapProducer", desc);
61 }
62 
64  auto hitMapHGCal = std::make_unique<DetIdRecHitMap>();
65  const auto& ee_hits = evt.get(hits_ee_token_);
66  const auto& fh_hits = evt.get(hits_fh_token_);
67  const auto& bh_hits = evt.get(hits_bh_token_);
68 
69  for (unsigned int i = 0; i < ee_hits.size(); ++i) {
70  hitMapHGCal->emplace(ee_hits[i].detid(), i);
71  }
72  auto size = ee_hits.size();
73 
74  for (unsigned int i = 0; i < fh_hits.size(); ++i) {
75  hitMapHGCal->emplace(fh_hits[i].detid(), i + size);
76  }
77  size += fh_hits.size();
78 
79  for (unsigned int i = 0; i < bh_hits.size(); ++i) {
80  hitMapHGCal->emplace(bh_hits[i].detid(), i + size);
81  }
82 
83  evt.put(std::move(hitMapHGCal), "hgcalRecHitMap");
84 
85  if (!hgcalOnly_) {
86  auto hitMapBarrel = std::make_unique<DetIdRecHitMap>();
87  const auto& eb_hits = evt.get(hits_eb_token_);
88  const auto& hb_hits = evt.get(hits_hb_token_);
89  const auto& ho_hits = evt.get(hits_ho_token_);
90  size = 0;
91 
92  for (unsigned int i = 0; i < eb_hits.size(); ++i) {
93  hitMapBarrel->emplace(eb_hits[i].detId(), i);
94  }
95  size += eb_hits.size();
96 
97  for (unsigned int i = 0; i < hb_hits.size(); ++i) {
98  hitMapBarrel->emplace(hb_hits[i].detId(), i + size);
99  }
100  size += hb_hits.size();
101 
102  for (unsigned int i = 0; i < ho_hits.size(); ++i) {
103  hitMapBarrel->emplace(ho_hits[i].detId(), i + size);
104  }
105  evt.put(std::move(hitMapBarrel), "barrelRecHitMap");
106  }
107 }
size
Write out results.
std::unordered_map< DetId, const unsigned int > DetIdRecHitMap
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
const edm::EDGetTokenT< reco::PFRecHitCollection > hits_ho_token_
bool get(ProductID const &oid, Handle< PROD > &result) const
Definition: Event.h:344
RecHitMapProducer(const edm::ParameterSet &)
const edm::EDGetTokenT< HGCRecHitCollection > hits_fh_token_
std::vector< PFRecHit > PFRecHitCollection
collection of PFRecHit objects
Definition: PFRecHitFwd.h:9
const edm::EDGetTokenT< reco::PFRecHitCollection > hits_hb_token_
const edm::EDGetTokenT< reco::PFRecHitCollection > hits_eb_token_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
const edm::EDGetTokenT< HGCRecHitCollection > hits_bh_token_
void add(std::string const &label, ParameterSetDescription const &psetDescription)
const edm::EDGetTokenT< HGCRecHitCollection > hits_ee_token_
fixed size matrix
HLT enums.
def move(src, dest)
Definition: eostools.py:511