CMS 3D CMS Logo

HitToLayerClusterAssociatorProducer.cc
Go to the documentation of this file.
1 // Author: Felice Pantaleo, felice.pantaleo@cern.ch 06/2024
2 
3 // user include files
15 
17  : LCCollectionToken_(consumes<std::vector<reco::CaloCluster>>(pset.getParameter<edm::InputTag>("layer_clusters"))),
18  hitMapToken_(consumes<std::unordered_map<DetId, unsigned int>>(pset.getParameter<edm::InputTag>("hitMap"))) {
19  auto hitsTags = pset.getParameter<std::vector<edm::InputTag>>("hits");
20  for (const auto &tag : hitsTags) {
21  hitsTokens_.push_back(consumes<HGCRecHitCollection>(tag));
22  }
23  produces<ticl::AssociationMap<ticl::mapWithFraction>>("hitToLayerClusterMap");
24 }
25 
27 
30  const edm::EventSetup &iSetup) const {
31  using namespace edm;
32 
35 
37  iEvent.getByToken(hitMapToken_, hitMap);
38 
39  MultiVectorManager<HGCRecHit> rechitManager;
40  for (const auto &token : hitsTokens_) {
41  Handle<HGCRecHitCollection> hitsHandle;
42  iEvent.getByToken(token, hitsHandle);
43  rechitManager.addVector(*hitsHandle);
44  }
45 
46  // Create association map
47  auto hitToLayerClusterMap = std::make_unique<ticl::AssociationMap<ticl::mapWithFraction>>(rechitManager.size());
48 
49  // Loop over layer clusters
50  for (unsigned int lcId = 0; lcId < layer_clusters->size(); ++lcId) {
51  const auto &layer_cluster = (*layer_clusters)[lcId];
52 
53  for (const auto &hitAndFraction : layer_cluster.hitsAndFractions()) {
54  auto hitMapIter = hitMap->find(hitAndFraction.first);
55  if (hitMapIter != hitMap->end()) {
56  unsigned int rechitIndex = hitMapIter->second;
57  float fraction = hitAndFraction.second;
58  hitToLayerClusterMap->insert(rechitIndex, lcId, fraction);
59  }
60  }
61  }
62  iEvent.put(std::move(hitToLayerClusterMap), "hitToLayerClusterMap");
63 }
64 
67  desc.add<edm::InputTag>("layer_clusters", edm::InputTag("hgcalMergeLayerClusters"));
68  desc.add<edm::InputTag>("hitMap", edm::InputTag("recHitMapProducer", "hgcalRecHitMap"));
69  desc.add<std::vector<edm::InputTag>>("hits",
70  {edm::InputTag("HGCalRecHit", "HGCEERecHits"),
71  edm::InputTag("HGCalRecHit", "HGCHEFRecHits"),
72  edm::InputTag("HGCalRecHit", "HGCHEBRecHits")});
73  descriptions.add("hitToLayerClusterAssociator", desc);
74 }
75 
76 // Define this as a plug-in
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::vector< edm::EDGetTokenT< HGCRecHitCollection > > hitsTokens_
Definition: DetId.h:17
void add(std::string const &label, ParameterSetDescription const &psetDescription)
void addVector(std::span< const T > vec)
fixed size matrix
edm::EDGetTokenT< std::vector< reco::CaloCluster > > LCCollectionToken_
HLT enums.
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
edm::EDGetTokenT< std::unordered_map< DetId, unsigned int > > hitMapToken_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
def move(src, dest)
Definition: eostools.py:511